CMS 3D CMS Logo

CosmicMuonTrajectoryBuilder.cc
Go to the documentation of this file.
12 
13 /* Collaborating Class Header */
37 
38 #include <algorithm>
39 
40 using namespace edm;
41 using namespace std;
42 
46  : theService(service) {
47  thePropagatorName = par.getParameter<string>("Propagator");
48 
49  bool enableDTMeasurement = par.getParameter<bool>("EnableDTMeasurement");
50  bool enableCSCMeasurement = par.getParameter<bool>("EnableCSCMeasurement");
51  bool enableRPCMeasurement = par.getParameter<bool>("EnableRPCMeasurement");
52 
53  // if(enableDTMeasurement)
54  InputTag DTRecSegmentLabel = par.getParameter<InputTag>("DTRecSegmentLabel");
55 
56  // if(enableCSCMeasurement)
57  InputTag CSCRecSegmentLabel = par.getParameter<InputTag>("CSCRecSegmentLabel");
58 
59  // if(enableRPCMeasurement)
60  InputTag RPCRecSegmentLabel = par.getParameter<InputTag>("RPCRecSegmentLabel");
61 
62  // if(enableGEMMeasurement)
63  // InputTag GEMRecSegmentLabel = par.getParameter<InputTag>("GEMRecSegmentLabel");
64 
65  theLayerMeasurements = new MuonDetLayerMeasurements(DTRecSegmentLabel,
66  CSCRecSegmentLabel,
67  RPCRecSegmentLabel,
68  edm::InputTag(),
69  edm::InputTag(),
70  iC,
71  enableDTMeasurement,
72  enableCSCMeasurement,
73  enableRPCMeasurement,
74  false,
75  false);
76 
77  ParameterSet muonUpdatorPSet = par.getParameter<ParameterSet>("MuonTrajectoryUpdatorParameters");
78 
79  theNavigation = nullptr; // new DirectMuonNavigation(theService->detLayerGeometry());
80  theUpdator = new MuonTrajectoryUpdator(muonUpdatorPSet, insideOut);
81 
83 
84  ParameterSet muonBackwardUpdatorPSet = par.getParameter<ParameterSet>("BackwardMuonTrajectoryUpdatorParameters");
85 
86  theBKUpdator = new MuonTrajectoryUpdator(muonBackwardUpdatorPSet, outsideIn);
87 
88  theTraversingMuonFlag = par.getParameter<bool>("BuildTraversingMuon");
89 
90  theStrict1LegFlag = par.getParameter<bool>("Strict1Leg");
91 
92  ParameterSet smootherPSet = par.getParameter<ParameterSet>("MuonSmootherParameters");
93 
94  theNavigationPSet = par.getParameter<ParameterSet>("MuonNavigationParameters");
95 
96  theSmoother = new CosmicMuonSmoother(smootherPSet, theService);
97 
98  theNTraversing = 0;
99  theNSuccess = 0;
100  theCacheId_DG = 0;
101  category_ = "Muon|RecoMuon|CosmicMuon|CosmicMuonTrajectoryBuilder";
102 }
103 
105  LogTrace(category_) << "CosmicMuonTrajectoryBuilder dtor called";
106  if (theUpdator)
107  delete theUpdator;
108  if (theBKUpdator)
109  delete theBKUpdator;
111  delete theLayerMeasurements;
112  if (theSmoother)
113  delete theSmoother;
114  if (theNavigation)
115  delete theNavigation;
117 
118  LogTrace(category_) << "CosmicMuonTrajectoryBuilder Traversing: " << theNSuccess << "/" << theNTraversing;
119 }
120 
123 
124  // DetLayer Geometry
125  unsigned long long newCacheId_DG = theService->eventSetup().get<MuonRecoGeometryRecord>().cacheIdentifier();
126  if (newCacheId_DG != theCacheId_DG) {
127  LogTrace(category_) << "Muon Reco Geometry changed!";
128  theCacheId_DG = newCacheId_DG;
129  if (theNavigation)
130  delete theNavigation;
132  }
133 }
134 
136  vector<Trajectory*> trajL = vector<Trajectory*>();
137 
139 
140  PTrajectoryStateOnDet ptsd1(seed.startingState());
141  DetId did(ptsd1.detId());
142  const BoundPlane& bp = theService->trackingGeometry()->idToDet(did)->surface();
143  TrajectoryStateOnSurface lastTsos =
144  trajectoryStateTransform::transientState(ptsd1, &bp, &*theService->magneticField());
145  LogTrace(category_) << "Seed: mom " << lastTsos.globalMomentum() << "pos: " << lastTsos.globalPosition();
146  LogTrace(category_) << "Seed: mom eta " << lastTsos.globalMomentum().eta()
147  << "pos eta: " << lastTsos.globalPosition().eta();
148 
149  bool beamhaloFlag = ((did.subdetId() == MuonSubdetId::CSC) && fabs(lastTsos.globalMomentum().eta()) > 4.0);
150 
151  vector<const DetLayer*> navLayers;
152 
153  if (did.subdetId() == MuonSubdetId::DT) {
154  //DT
155  navLayers = navigation()->compatibleLayers(*(lastTsos.freeState()), alongMomentum);
156  } else if (beamhaloFlag || (theTraversingMuonFlag && theStrict1LegFlag)) {
157  //CSC
158  navLayers = navigation()->compatibleEndcapLayers(*(lastTsos.freeState()), alongMomentum);
159  } else {
160  navLayers = navigation()->compatibleLayers(*(lastTsos.freeState()), alongMomentum);
161  }
162 
163  LogTrace(category_) << "found " << navLayers.size() << " compatible DetLayers for the Seed";
164 
165  if (navLayers.empty())
166  return trajL;
167 
168  vector<DetWithState> detsWithStates;
169  LogTrace(category_) << "Compatible layers: ";
170  for (vector<const DetLayer*>::const_iterator layer = navLayers.begin(); layer != navLayers.end(); layer++) {
171  LogTrace(category_) << debug.dumpMuonId((*layer)->basicComponents().front()->geographicalId())
172  << debug.dumpLayer(*layer);
173  }
174 
175  detsWithStates = navLayers.front()->compatibleDets(lastTsos, *propagator(), *(updator()->estimator()));
176  LogTrace(category_) << "Number of compatible dets: " << detsWithStates.size() << endl;
177 
178  if (!detsWithStates.empty()) {
179  // get the updated TSOS
180  if (detsWithStates.front().second.isValid()) {
181  LogTrace(category_) << "New starting TSOS is on det: " << endl;
182  LogTrace(category_) << debug.dumpMuonId(detsWithStates.front().first->geographicalId())
183  << debug.dumpLayer(navLayers.front());
184  lastTsos = detsWithStates.front().second;
185  LogTrace(category_) << "Seed after extrapolation: mom " << lastTsos.globalMomentum()
186  << "pos: " << lastTsos.globalPosition();
187  }
188  }
189  detsWithStates.clear();
190  if (!lastTsos.isValid())
191  return trajL;
192 
193  TrajectoryStateOnSurface secondLast = lastTsos;
194 
195  lastTsos.rescaleError(10.0);
196 
197  Trajectory* theTraj = new Trajectory(seed, alongMomentum);
198 
199  navLayers.clear();
200 
201  if (fabs(lastTsos.globalMomentum().eta()) < 1.0) {
202  //DT
203  navLayers = navigation()->compatibleLayers(*(lastTsos.freeState()), alongMomentum);
204  } else if (beamhaloFlag || (theTraversingMuonFlag && theStrict1LegFlag)) {
205  //CSC
206  navLayers = navigation()->compatibleEndcapLayers(*(lastTsos.freeState()), alongMomentum);
207  } else {
208  navLayers = navigation()->compatibleLayers(*(lastTsos.freeState()), alongMomentum);
209  }
210 
211  int DTChamberUsedBack = 0;
212  int CSCChamberUsedBack = 0;
213  int RPCChamberUsedBack = 0;
214  int TotalChamberUsedBack = 0;
216  vector<TrajectoryMeasurement> measL;
217 
218  LogTrace(category_) << "Begin forward fit " << navLayers.size();
219 
220  for (vector<const DetLayer*>::const_iterator rnxtlayer = navLayers.begin(); rnxtlayer != navLayers.end();
221  ++rnxtlayer) {
222  LogTrace(category_) << "new layer ";
223  measL.clear();
224  LogTrace(category_) << debug.dumpMuonId((*rnxtlayer)->basicComponents().front()->geographicalId())
225  << debug.dumpLayer(*rnxtlayer);
226  LogTrace(category_) << "from lastTsos " << lastTsos.globalMomentum() << " at " << lastTsos.globalPosition();
227 
228  measL = findBestMeasurements(*rnxtlayer, lastTsos, propagator(), (updator()->estimator()));
229 
230  if (measL.empty() && (fabs(theService->magneticField()->inTesla(GlobalPoint(0, 0, 0)).z()) < 0.01) &&
231  (theService->propagator("StraightLinePropagator").isValid())) {
232  LogTrace(category_) << "try straight line propagator ";
233  measL = findBestMeasurements(
234  *rnxtlayer, lastTsos, &*theService->propagator("StraightLinePropagator"), (updator()->estimator()));
235  }
236  if (measL.empty())
237  continue;
238 
239  for (vector<TrajectoryMeasurement>::const_iterator theMeas = measL.begin(); theMeas != measL.end(); ++theMeas) {
240  pair<bool, TrajectoryStateOnSurface> result = updator()->update((&*theMeas), *theTraj, propagator());
241 
242  if (result.first) {
243  LogTrace(category_) << "update ok ";
245  (*rnxtlayer), DTChamberUsedBack, CSCChamberUsedBack, RPCChamberUsedBack, TotalChamberUsedBack);
246  secondLast = lastTsos;
247  if ((!theTraj->empty()) && result.second.isValid()) {
248  lastTsos = result.second;
249  LogTrace(category_) << "get new lastTsos here " << lastTsos.globalMomentum() << " at "
250  << lastTsos.globalPosition();
251  } else if ((theMeas)->predictedState().isValid())
252  lastTsos = (theMeas)->predictedState();
253  }
254  }
255  }
256  measL.clear();
257  while (!theTraj->empty()) {
258  theTraj->pop();
259  }
260 
261  if (!theTraj->isValid() || TotalChamberUsedBack < 2 || (DTChamberUsedBack + CSCChamberUsedBack) == 0 ||
262  !lastTsos.isValid()) {
263  delete theTraj;
264  return trajL;
265  }
266  delete theTraj;
267 
268  // if got good trajectory, then do backward refitting
269  DTChamberUsedBack = 0;
270  CSCChamberUsedBack = 0;
271  RPCChamberUsedBack = 0;
272  TotalChamberUsedBack = 0;
273 
274  Trajectory myTraj(seed, oppositeToMomentum);
275 
276  // set starting navigation direction for MuonTrajectoryUpdator
277 
278  GlobalPoint lastPos = lastTsos.globalPosition();
279  GlobalPoint secondLastPos = secondLast.globalPosition();
280  GlobalVector momDir = secondLastPos - lastPos;
281 
282  if (lastPos.basicVector().dot(momDir.basicVector()) > 0) {
283  // LogTrace("CosmicMuonTrajectoryBuilder")<<"Fit direction changed to insideOut";
285  } else
287 
288  if (fabs(lastTsos.globalMomentum().eta()) < 1.0) {
289  //DT
290  navLayers = navigation()->compatibleLayers(*(lastTsos.freeState()), oppositeToMomentum);
291  } else if (beamhaloFlag || (theTraversingMuonFlag && theStrict1LegFlag)) {
292  //CSC
293  std::reverse(navLayers.begin(), navLayers.end());
294  } else {
295  navLayers = navigation()->compatibleLayers(*(lastTsos.freeState()), oppositeToMomentum);
296  }
297 
298  LogTrace(category_) << "Begin backward refitting, with " << navLayers.size() << " layers" << endl;
299 
300  for (vector<const DetLayer*>::const_iterator rnxtlayer = navLayers.begin(); rnxtlayer != navLayers.end();
301  ++rnxtlayer) {
302  measL.clear();
303 
304  measL = findBestMeasurements(*rnxtlayer, lastTsos, propagator(), (backwardUpdator()->estimator()));
305 
306  if (measL.empty()) {
308  for (MuonRecHitContainer::const_iterator ihit = tmpHits.begin(); ihit != tmpHits.end(); ++ihit) {
309  allUnusedHits.push_back(*ihit);
310  }
311  continue;
312  }
313 
314  for (vector<TrajectoryMeasurement>::const_iterator theMeas = measL.begin(); theMeas != measL.end(); ++theMeas) {
315  // if the part change, we need to reconsider the fit direction
316  if (rnxtlayer != navLayers.begin()) {
317  vector<const DetLayer*>::const_iterator lastlayer = rnxtlayer;
318  lastlayer--;
319 
320  if ((*rnxtlayer)->location() != (*lastlayer)->location()) {
321  lastPos = lastTsos.globalPosition();
322  GlobalPoint thisPos = (theMeas)->predictedState().globalPosition();
323  GlobalVector momDir = thisPos - lastPos;
324  // LogTrace("CosmicMuonTrajectoryBuilder")<<"momDir "<<momDir;
325 
326  if (momDir.mag() > 0.01) { //if lastTsos is on the surface, no need
327  if (thisPos.basicVector().dot(momDir.basicVector()) > 0) {
329  } else
331  }
332  }
333  if (((*lastlayer)->location() == GeomDetEnumerators::endcap) &&
334  ((*rnxtlayer)->location() == GeomDetEnumerators::endcap) &&
335  (lastTsos.globalPosition().z() * (theMeas)->predictedState().globalPosition().z() < 0)) {
337  }
338  }
339 
340  // if (theBKUpdator->fitDirection() == insideOut)
341  // LogTrace("CosmicMuonTrajectoryBuilder")<<"Fit direction insideOut";
342  // else LogTrace("CosmicMuonTrajectoryBuilder")<<"Fit direction outsideIn";
343  pair<bool, TrajectoryStateOnSurface> bkresult = backwardUpdator()->update((&*theMeas), myTraj, propagator());
344 
345  if (bkresult.first) {
347  (*rnxtlayer), DTChamberUsedBack, CSCChamberUsedBack, RPCChamberUsedBack, TotalChamberUsedBack);
348 
349  if (theTraversingMuonFlag) {
350  MuonRecHitContainer tmpUnusedHits = unusedHits(*rnxtlayer, *theMeas);
351  allUnusedHits.insert(allUnusedHits.end(), tmpUnusedHits.begin(), tmpUnusedHits.end());
352  }
353  if ((!myTraj.empty()) && bkresult.second.isValid())
354  lastTsos = bkresult.second;
355  else if ((theMeas)->predictedState().isValid())
356  lastTsos = (theMeas)->predictedState();
357  }
358  }
359  }
360 
361  for (vector<Trajectory*>::iterator t = trajL.begin(); t != trajL.end(); ++t)
362  delete *t;
363 
364  trajL.clear();
365 
366  if ((!myTraj.isValid()) || (myTraj.empty()) || ((selfDuplicate(myTraj))) || TotalChamberUsedBack < 2 ||
367  (DTChamberUsedBack + CSCChamberUsedBack) < 1) {
368  return trajL;
369  }
370 
371  if (theTraversingMuonFlag && (allUnusedHits.size() >= 2)) {
372  // LogTrace(category_)<<utilities()->print(allUnusedHits);
373  LogTrace(category_) << "Building trajectory in second hemisphere...";
374  buildSecondHalf(myTraj);
375  // check if traversing trajectory has hits in both hemispheres
376 
377  if (theStrict1LegFlag && !utilities()->isTraversing(myTraj)) {
378  trajL.clear();
379  return trajL;
380  }
382  trajL.clear();
383  return trajL;
384  }
385 
386  LogTrace(category_) << " traj ok ";
387 
388  // getDirectionByTime(myTraj);
389  if (beamhaloFlag)
390  estimateDirection(myTraj);
391  if (myTraj.empty())
392  return trajL;
393 
394  // try to smooth it
395  vector<Trajectory> smoothed = theSmoother->trajectories(myTraj);
396 
397  if (!smoothed.empty() && smoothed.front().foundHits() > 3) {
398  LogTrace(category_) << " Smoothed successfully.";
399  myTraj = smoothed.front();
400  } else {
401  LogTrace(category_) << " Smooth failed.";
402  }
403 
404  LogTrace(category_) << "first " << myTraj.firstMeasurement().updatedState() << "\n last "
405  << myTraj.lastMeasurement().updatedState();
406  if (myTraj.direction() == alongMomentum)
407  LogTrace(category_) << "alongMomentum";
408  else if (myTraj.direction() == oppositeToMomentum)
409  LogTrace(category_) << "oppositeMomentum";
410  else
411  LogTrace(category_) << "anyDirection";
412 
413  if (!beamhaloFlag) {
414  if (myTraj.lastMeasurement().updatedState().globalMomentum().y() > 0) {
415  LogTrace(category_) << "flip trajectory ";
416  flipTrajectory(myTraj);
417  }
418 
419  if ((myTraj.direction() == alongMomentum && (myTraj.firstMeasurement().updatedState().globalPosition().y() <
420  myTraj.lastMeasurement().updatedState().globalPosition().y())) ||
422  myTraj.lastMeasurement().updatedState().globalPosition().y()))) {
423  LogTrace(category_) << "reverse propagation direction";
425  }
426  }
427  // getDirectionByTime(myTraj);
428  if (!myTraj.isValid())
429  return trajL;
430 
431  // check direction agree with position!
433  GlobalVector dirFromPos = myTraj.measurements().back().recHit()->globalPosition() -
434  myTraj.measurements().front().recHit()->globalPosition();
435 
436  if (theStrict1LegFlag && !utilities()->isTraversing(myTraj)) {
437  trajL.clear();
438  return trajL;
439  }
440 
441  LogTrace(category_) << "last hit " << myTraj.measurements().back().recHit()->globalPosition() << endl;
442  LogTrace(category_) << "first hit " << myTraj.measurements().front().recHit()->globalPosition() << endl;
443 
444  LogTrace(category_) << "last tsos " << myTraj.measurements().back().updatedState().globalPosition() << " mom "
445  << myTraj.measurements().back().updatedState().globalMomentum() << endl;
446  LogTrace(category_) << "first tsos " << myTraj.measurements().front().updatedState().globalPosition() << " mom "
447  << myTraj.measurements().front().updatedState().globalMomentum() << endl;
448 
449  PropagationDirection propDir =
450  (dirFromPos.basicVector().dot(myTraj.firstMeasurement().updatedState().globalMomentum().basicVector()) > 0)
451  ? alongMomentum
453  LogTrace(category_) << " dir " << dir << " propDir " << propDir << endl;
454 
455  LogTrace(category_) << "chi2 " << myTraj.chiSquared() << endl;
456 
457  if (dir != propDir) {
458  LogTrace(category_) << "reverse propagation direction ";
460  }
461  if (myTraj.empty())
462  return trajL;
463 
464  trajL.push_back(new Trajectory(myTraj));
465  navLayers.clear();
466  return trajL;
467 }
468 
469 //
470 //
471 //
473  const DetLayer* layer, const TrajectoryMeasurement& meas) const {
476  for (MuonRecHitContainer::const_iterator ihit = tmpHits.begin(); ihit != tmpHits.end(); ++ihit) {
477  if ((*ihit)->geographicalId() != meas.recHit()->geographicalId()) {
478  result.push_back(*ihit);
479  LogTrace(category_) << "Unused hit: " << (*ihit)->globalPosition() << endl;
480  }
481  }
482 
483  return result;
484 }
485 
486 //
487 // continue to build a trajectory starting from given trajectory state
488 //
490  const NavigationDirection& startingDir,
491  Trajectory& traj) {
492  if (!ts.isValid())
493  return;
494 
495  const FreeTrajectoryState* fts = ts.freeState();
496  if (!fts)
497  return;
498 
499  vector<const DetLayer*> navLayers;
500 
501  if (fabs(fts->momentum().basicVector().eta()) < 1.0) {
502  //DT
503  if (fts->position().basicVector().dot(fts->momentum().basicVector()) > 0) {
504  navLayers = navigation()->compatibleLayers((*fts), alongMomentum);
505  } else {
506  navLayers = navigation()->compatibleLayers((*fts), oppositeToMomentum);
507  }
508 
510  //CSC
511  if (fts->position().basicVector().dot(fts->momentum().basicVector()) > 0) {
512  navLayers = navigation()->compatibleEndcapLayers((*fts), alongMomentum);
513  } else {
514  navLayers = navigation()->compatibleEndcapLayers((*fts), oppositeToMomentum);
515  }
516  } else {
517  if (fts->position().basicVector().dot(fts->momentum().basicVector()) > 0) {
518  navLayers = navigation()->compatibleLayers((*fts), alongMomentum);
519  } else {
520  navLayers = navigation()->compatibleLayers((*fts), oppositeToMomentum);
521  }
522  }
523 
524  if (navLayers.empty())
525  return;
526 
527  theBKUpdator->setFitDirection(startingDir);
528 
529  int DTChamberUsedBack = 0;
530  int CSCChamberUsedBack = 0;
531  int RPCChamberUsedBack = 0;
532  int TotalChamberUsedBack = 0;
533 
536  ? propagatorAlong()->propagate((*fts), navLayers.front()->surface())
537  : propagatorOpposite()->propagate((*fts), navLayers.front()->surface());
538 
539  if (!lastTsos.isValid()) {
540  LogTrace(category_) << "propagation failed from fts to inner cylinder";
541  return;
542  }
543  LogTrace(category_) << "tsos " << lastTsos.globalPosition();
544  lastTsos.rescaleError(10.);
545  vector<TrajectoryMeasurement> measL;
546  for (vector<const DetLayer*>::const_iterator rnxtlayer = navLayers.begin(); rnxtlayer != navLayers.end();
547  ++rnxtlayer) {
548  measL.clear();
549  measL = findBestMeasurements(*rnxtlayer, lastTsos, propagator(), (backwardUpdator()->estimator()));
550 
551  if (measL.empty())
552  continue;
553 
554  for (vector<TrajectoryMeasurement>::const_iterator theMeas = measL.begin(); theMeas != measL.end(); ++theMeas) {
555  pair<bool, TrajectoryStateOnSurface> bkresult = backwardUpdator()->update((&*theMeas), traj, propagator());
556  if (bkresult.first) {
557  LogTrace(category_) << "update ok : " << (theMeas)->recHit()->globalPosition();
558 
560  (*rnxtlayer), DTChamberUsedBack, CSCChamberUsedBack, RPCChamberUsedBack, TotalChamberUsedBack);
561 
562  if ((!traj.empty()) && bkresult.second.isValid())
563  lastTsos = bkresult.second;
564  else if ((theMeas)->predictedState().isValid())
565  lastTsos = (theMeas)->predictedState();
566  }
567  }
568  }
569  navLayers.clear();
570  updator()->makeFirstTime();
572 
573  measL.clear();
574 
575  return;
576 }
577 
578 //
579 // build trajectory in second hemisphere with pattern
580 // recognition starting from an intermediate state
581 //
583  if ((traj.firstMeasurement().recHit()->globalPosition().perp() <
584  traj.lastMeasurement().recHit()->globalPosition().perp())) {
585  LogTrace(category_) << "inside-out: reverseTrajectory";
586  reverseTrajectory(traj);
587  }
588  if (traj.empty())
589  return;
591  if (!tsos.isValid())
592  tsos = traj.lastMeasurement().predictedState();
593  LogTrace(category_) << "last tsos on traj: pos: " << tsos.globalPosition() << " mom: " << tsos.globalMomentum();
594  if (!tsos.isValid()) {
595  LogTrace(category_) << "last tsos on traj invalid";
596  return;
597  }
598 
599  build(intermediateState(tsos), insideOut, traj);
600 
601  return;
602 }
603 
604 //
605 //
606 //
608  PerpendicularBoundPlaneBuilder planeBuilder;
609  GlobalPoint pos(0.0, 0.0, 0.0);
610  BoundPlane* SteppingPlane = planeBuilder(pos, tsos.globalDirection());
611 
612  TrajectoryStateOnSurface predTsos = propagator()->propagate(tsos, *SteppingPlane);
613  if (predTsos.isValid())
614  LogTrace(category_) << "intermediateState: a intermediate state: pos: " << predTsos.globalPosition()
615  << "mom: " << predTsos.globalMomentum();
616 
617  return predTsos;
618 }
619 
620 //
621 //
622 //
624  if (hits.size() < 2)
625  return;
626 
628  vector<bool> keep(hits.size(), true);
629  int i(0);
630  int j(0);
631 
632  for (MuonRecHitContainer::const_iterator ihit = hits.begin(); ihit != hits.end(); ++ihit) {
633  if (!keep[i]) {
634  i++;
635  continue;
636  };
637  j = i + 1;
638  for (MuonRecHitContainer::const_iterator ihit2 = ihit + 1; ihit2 != hits.end(); ++ihit2) {
639  if (!keep[j]) {
640  j++;
641  continue;
642  }
643  if ((*ihit)->geographicalId() == (*ihit2)->geographicalId()) {
644  if ((*ihit)->dimension() > (*ihit2)->dimension()) {
645  keep[j] = false;
646  } else if ((*ihit)->dimension() < (*ihit2)->dimension()) {
647  keep[i] = false;
648  } else {
649  if ((*ihit)->transientHits().size() > (*ihit2)->transientHits().size()) {
650  keep[j] = false;
651  } else if ((*ihit)->transientHits().size() < (*ihit2)->transientHits().size()) {
652  keep[i] = false;
653  } else if ((*ihit)->degreesOfFreedom() != 0 && (*ihit2)->degreesOfFreedom() != 0) {
654  if (((*ihit)->chi2() / (*ihit)->degreesOfFreedom()) > ((*ihit2)->chi2() / (*ihit)->degreesOfFreedom()))
655  keep[i] = false;
656  else
657  keep[j] = false;
658  }
659  }
660  } // if same geomid
661  j++;
662  }
663  i++;
664  }
665 
666  i = 0;
667  for (MuonRecHitContainer::const_iterator ihit = hits.begin(); ihit != hits.end(); ++ihit) {
668  if (keep[i])
669  tmp.push_back(*ihit);
670  i++;
671  }
672 
673  hits.clear();
674  hits.swap(tmp);
675 
676  return;
677 }
678 
679 //
680 //
681 //
684 
685  if (traj.empty())
686  return true;
687 
688  bool result = false;
689  for (ConstRecHitContainer::const_iterator ir = hits.begin(); ir != hits.end(); ir++) {
690  if (!(*ir)->isValid())
691  continue;
692  for (ConstRecHitContainer::const_iterator ir2 = ir + 1; ir2 != hits.end(); ir2++) {
693  if (!(*ir2)->isValid())
694  continue;
695  if ((*ir) == (*ir2))
696  result = true;
697  }
698  }
699 
700  return result;
701 }
702 
703 //
704 // reverse a trajectory without refitting
705 // this can only be used for cosmic muons that come from above
706 //
708  PropagationDirection newDir =
709  (traj.firstMeasurement().recHit()->globalPosition().y() < traj.lastMeasurement().recHit()->globalPosition().y())
711  : alongMomentum;
712  Trajectory newTraj(traj.seed(), newDir);
713 
714  /* does not work in gcc4.8?)
715  std::vector<TrajectoryMeasurement> & meas = traj.measurements();
716  for (auto itm = meas.rbegin(); itm != meas.rend(); ++itm ) {
717  newTraj.push(std::move(*itm));
718  }
719  traj = std::move(newTraj);
720  */
721 
722  std::vector<TrajectoryMeasurement> const& meas = traj.measurements();
723  for (auto itm = meas.rbegin(); itm != meas.rend(); ++itm) {
724  newTraj.push(*itm);
725  }
726  traj = newTraj;
727 }
728 
729 //
730 // reverse a trajectory momentum direction and then refit
731 //
734  if (!lastTSOS.isValid()) {
735  LogTrace(category_) << "Error: last TrajectoryState invalid.";
736  }
738  std::reverse(hits.begin(), hits.end());
739 
740  LogTrace(category_) << "last tsos before flipping " << lastTSOS;
741  utilities()->reverseDirection(lastTSOS, &*theService->magneticField());
742  LogTrace(category_) << "last tsos after flipping " << lastTSOS;
743 
744  vector<Trajectory> refittedback = theSmoother->fit(traj.seed(), hits, lastTSOS);
745  if (refittedback.empty()) {
746  LogTrace(category_) << "flipTrajectory fail. " << endl;
747  return;
748  }
749  LogTrace(category_) << "flipTrajectory: first " << refittedback.front().firstMeasurement().updatedState()
750  << "\nflipTrajectory: last " << refittedback.front().lastMeasurement().updatedState();
751 
752  traj = refittedback.front();
753 
754  return;
755 }
756 
757 //
758 //
759 //
761  if (traj.direction() == anyDirection)
762  return;
764  Trajectory newTraj(traj.seed(), newDir);
765  const std::vector<TrajectoryMeasurement>& meas = traj.measurements();
766 
767  for (std::vector<TrajectoryMeasurement>::const_iterator itm = meas.begin(); itm != meas.end(); ++itm) {
768  newTraj.push(*itm);
769  }
770 
771  while (!traj.empty()) {
772  traj.pop();
773  }
774 
775  traj = newTraj;
776 }
777 
778 //
779 // guess the direction by normalized chi2
780 //
783 
785 
787 
788  if (!firstTSOS.isValid() || !lastTSOS.isValid())
789  return;
790 
791  LogTrace(category_) << "Two ends of the traj " << firstTSOS.globalPosition() << ", " << lastTSOS.globalPosition();
792 
793  LogTrace(category_) << "Their mom: " << firstTSOS.globalMomentum() << ", " << lastTSOS.globalMomentum();
794 
795  LogTrace(category_) << "Their mom eta: " << firstTSOS.globalMomentum().eta() << ", "
796  << lastTSOS.globalMomentum().eta();
797 
798  // momentum eta can be used to estimate direction
799  // the beam-halo muon seems enter with a larger |eta|
800 
801  if (fabs(firstTSOS.globalMomentum().eta()) > fabs(lastTSOS.globalMomentum().eta())) {
802  vector<Trajectory> refitted = theSmoother->trajectories(traj.seed(), hits, firstTSOS);
803  if (!refitted.empty())
804  traj = refitted.front();
805 
806  } else {
807  std::reverse(hits.begin(), hits.end());
808  utilities()->reverseDirection(lastTSOS, &*theService->magneticField());
809  vector<Trajectory> refittedback = theSmoother->trajectories(traj.seed(), hits, lastTSOS);
810  if (!refittedback.empty())
811  traj = refittedback.front();
812  }
813 
814  return;
815 }
816 
817 //
818 // get direction from timing information of rechits and segments
819 //
822  LogTrace(category_) << "getDirectionByTime" << endl;
823  for (TransientTrackingRecHit::ConstRecHitContainer::const_iterator ir = hits.begin(); ir != hits.end(); ir++) {
824  if (!(*ir)->isValid()) {
825  LogTrace(category_) << "invalid RecHit" << endl;
826  continue;
827  }
828 
829  const GlobalPoint& pos = (*ir)->globalPosition();
830  LogTrace(category_) << "pos" << pos << "radius " << pos.perp() << " dim " << (*ir)->dimension() << " det "
831  << (*ir)->det()->geographicalId().det() << " sub det " << (*ir)->det()->subDetector() << endl;
832 
833  if ((*ir)->det()->geographicalId().det() == 2 && (*ir)->det()->subDetector() == 6) {
834  // const CSCRecHit2D* iCSC = dynamic_cast<const CSCRecHit2D*>(&**ir);
835  // if (iCSC) LogTrace(category_)<<"csc from cast tpeak "<<iCSC->tpeak();
836  CSCRecHit2DCollection::range thisrange = cschits_->get(CSCDetId((*ir)->geographicalId()));
837  for (CSCRecHit2DCollection::const_iterator rechit = thisrange.first; rechit != thisrange.second; ++rechit) {
838  if ((*rechit).isValid())
839  LogTrace(category_) << "csc from collection tpeak " << (*rechit).tpeak();
840  }
841  }
842  if ((*ir)->det()->geographicalId().det() == 2 && (*ir)->det()->subDetector() == 7) {
843  // const DTRecHit1D* iDT = dynamic_cast<const DTRecHit1D*>(&**ir);
844  // if (iDT) LogTrace(category_)<<"dt digitime "<<iDT->digiTime();
845  DTRecHitCollection::range thisrange = dthits_->get(DTLayerId((*ir)->geographicalId()));
846  for (DTRecHitCollection::const_iterator rechit = thisrange.first; rechit != thisrange.second; ++rechit) {
847  if ((*rechit).isValid())
848  LogTrace(category_) << "dt from collection digitime " << (*rechit).digiTime();
849  }
850  }
851  }
852 
853  return;
854 }
855 
856 //
857 //
858 //
859 std::vector<TrajectoryMeasurement> CosmicMuonTrajectoryBuilder::findBestMeasurements(
860  const DetLayer* layer,
861  const TrajectoryStateOnSurface& tsos,
862  const Propagator* propagator,
864  std::vector<TrajectoryMeasurement> result;
865  std::vector<TrajectoryMeasurement> measurements;
866 
867  if (layer->hasGroups()) {
868  std::vector<TrajectoryMeasurementGroup> measurementGroups =
869  theLayerMeasurements->groupedMeasurements(layer, tsos, *propagator, *estimator);
870 
871  for (std::vector<TrajectoryMeasurementGroup>::const_iterator tmGroupItr = measurementGroups.begin();
872  tmGroupItr != measurementGroups.end();
873  ++tmGroupItr) {
874  measurements = tmGroupItr->measurements();
875  const TrajectoryMeasurement* bestMeasurement =
876  theBestMeasurementFinder->findBestMeasurement(measurements, propagator);
877 
878  if (bestMeasurement)
879  result.push_back(*bestMeasurement);
880  }
881  } else {
882  measurements = theLayerMeasurements->measurements(layer, tsos, *propagator, *estimator);
883  const TrajectoryMeasurement* bestMeasurement =
884  theBestMeasurementFinder->findBestMeasurement(measurements, propagator);
885 
886  if (bestMeasurement)
887  result.push_back(*bestMeasurement);
888  }
889  measurements.clear();
890 
891  return result;
892 }
893 
894 //
895 //
896 //
898  const DetLayer* layer, int& dtChambers, int& cscChambers, int& rpcChambers, int& totalChambers) {
899  if (layer->subDetector() == GeomDetEnumerators::DT)
900  dtChambers++;
901  else if (layer->subDetector() == GeomDetEnumerators::CSC)
902  cscChambers++;
903  else if (layer->subDetector() == GeomDetEnumerators::RPCBarrel ||
905  rpcChambers++;
906  totalChambers++;
907 }
908 
909 //
910 //
911 //
913  if ((dtseg == nullptr) || (!dtseg->hasPhi()))
914  return 0;
915  // timing information
916  double result = 0;
917  if (dtseg->phiSegment() == nullptr)
918  return 0;
919  int phiHits = dtseg->phiSegment()->specificRecHits().size();
920  LogTrace(category_) << "phiHits " << phiHits;
921  if (phiHits > 5) {
922  if (dtseg->phiSegment()->ist0Valid())
923  result = dtseg->phiSegment()->t0();
924  if (dtseg->phiSegment()->ist0Valid()) {
925  LogTrace(category_) << " Phi t0: " << dtseg->phiSegment()->t0() << " hits: " << phiHits;
926  } else {
927  LogTrace(category_) << " Phi t0 is invalid: " << dtseg->phiSegment()->t0() << " hits: " << phiHits;
928  }
929  }
930 
931  return result;
932 }
933 
934 //
935 //
936 //
938  const DTRecSegment4D* dtseg2) const {
939  LogTrace(category_) << "comparing dtseg: " << dtseg1 << " " << dtseg2 << endl;
940  if (dtseg1 == dtseg2 || t0(dtseg1) == t0(dtseg2))
941  return anyDirection;
942 
944 
945  return result;
946 }
bool empty() const
True if trajectory has no measurements.
Definition: Trajectory.h:233
void reverseTrajectoryPropagationDirection(Trajectory &) const
reverse the propagation direction of a trajectory
T getParameter(std::string const &) const
const Propagator * propagator() const
double t0(const DTRecSegment4D *deseg) const
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
TrajectoryStateOnSurface const & predictedState() const
const MeasurementEstimator * estimator() const
accasso at the propagator
T perp() const
Definition: PV3DBase.h:69
ConstRecHitPointer const & recHit() const
TrajectorySeed const & seed() const
Access to the seed used to reconstruct the Trajectory.
Definition: Trajectory.h:263
std::vector< const DetLayer * > compatibleEndcapLayers(const FreeTrajectoryState &fts, PropagationDirection timeDirection) const
void estimateDirection(Trajectory &) const
check the direction of trajectory by checking eta spread
void makeFirstTime()
reset the theFirstTSOSFlag
TrajectoryStateOnSurface intermediateState(const TrajectoryStateOnSurface &) const
edm::Handle< CSCRecHit2DCollection > cschits_
virtual SubDetector subDetector() const =0
The type of detector (PixelBarrel, PixelEndcap, TIB, TOB, TID, TEC, CSC, DT, RPCBarrel, RPCEndcap)
const DTChamberRecSegment2D * phiSegment() const
The superPhi segment: 0 if no phi projection available.
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:60
MuonTrajectoryUpdator * backwardUpdator() const
GlobalPoint globalPosition() const
ConstRecHitContainer recHits() const
Definition: Trajectory.h:186
MuonCandidate::TrajectoryContainer TrajectoryContainer
edm::Handle< DTRecHitCollection > dthits_
const Propagator * propagatorOpposite() const
void setEvent(const edm::Event &) override
pass the Event to the algo at each event
PropagationDirection
MuonDetLayerMeasurements * theLayerMeasurements
const CosmicMuonUtilities * utilities() const
std::vector< TrajectoryMeasurementGroup > groupedMeasurements(const DetLayer *layer, const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est, const edm::Event &iEvent)
std::vector< Trajectory * > trajectories(const TrajectorySeed &) override
build trajectories from seed
std::vector< const DetLayer * > compatibleLayers(const FreeTrajectoryState &fts, PropagationDirection timeDirection) const
PropagationDirection const & direction() const
Definition: Trajectory.cc:133
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
DataContainer const & measurements() const
Definition: Trajectory.h:178
const int keep
DirectMuonNavigation * navigation() const
T mag() const
Definition: PV3DBase.h:64
void flipTrajectory(Trajectory &) const
flip a trajectory with refit (the momentum direction is opposite)
MuonTrajectoryUpdator * updator() const
MuonRecHitContainer recHits(const DetLayer *layer, const edm::Event &iEvent)
returns the rechits which are on the layer
void reverseTrajectory(Trajectory &) const
reverse a trajectory without refit (out the measurements order changed)
MuonBestMeasurementFinder * theBestMeasurementFinder
void setFitDirection(NavigationDirection fitDirection)
set fit direction
TrajectoryMeasurement const & lastMeasurement() const
Definition: Trajectory.h:150
FreeTrajectoryState const * freeState(bool withErrors=true) const
T z() const
Definition: PV3DBase.h:61
void incrementChamberCounters(const DetLayer *layer, int &dtChambers, int &cscChambers, int &rpcChambers, int &totalChambers)
bool hasPhi() const
Does it have the Phi projection?
~CosmicMuonTrajectoryBuilder() override
Destructor.
MeasurementContainer measurements(const DetLayer *layer, const GeomDet *det, const TrajectoryStateOnSurface &stateOnDet, const MeasurementEstimator &est, const edm::Event &iEvent)
GlobalVector momentum() const
#define LogTrace(id)
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
bool ist0Valid() const
std::vector< ConstRecHitPointer > ConstRecHitContainer
void getDirectionByTime(Trajectory &) const
check the direction of trajectory by checking the timing
Definition: DetId.h:17
GlobalPoint position() const
PropagationDirection checkDirectionByT0(const DTRecSegment4D *, const DTRecSegment4D *) const
bool isValid() const
Definition: Trajectory.h:257
TrajectoryMeasurement const & firstMeasurement() const
Definition: Trajectory.h:166
PTrajectoryStateOnDet const & startingState() const
#define debug
Definition: HDRShower.cc:19
TrajectoryStateOnSurface transientState(const PTrajectoryStateOnDet &ts, const Surface *surface, const MagneticField *field)
void pop()
Definition: Trajectory.cc:30
std::vector< Trajectory > fit(const Trajectory &) const
TrajectoryStateOnSurface propagate(STA const &state, SUR const &surface) const
Definition: Propagator.h:50
float chiSquared() const
Definition: Trajectory.h:241
T eta() const
Definition: PV3DBase.h:73
void selectHits(MuonTransientTrackingRecHit::MuonRecHitContainer &) const
void setEvent(const edm::Event &)
set event
HLT enums.
GlobalVector globalMomentum() const
CosmicMuonTrajectoryBuilder(const edm::ParameterSet &, const MuonServiceProxy *service, edm::ConsumesCollector &iC)
Constructor.
const Propagator * propagatorAlong() const
TrajectoryMeasurement * findBestMeasurement(std::vector< TrajectoryMeasurement > &measC, const Propagator *propagator)
return the Tm with the best chi2: no cut applied.
MuonTransientTrackingRecHit::MuonRecHitContainer MuonRecHitContainer
std::vector< TrajectoryMeasurement > findBestMeasurements(const DetLayer *, const TrajectoryStateOnSurface &, const Propagator *, const MeasurementEstimator *)
virtual std::pair< bool, TrajectoryStateOnSurface > update(const TrajectoryMeasurement *measurement, Trajectory &trajectory, const Propagator *propagator)
update the Trajectory with the TrajectoryMeasurement
bool selfDuplicate(const Trajectory &) const
check if the trajectory iterates the same hit more than once
TrajectoryStateOnSurface const & updatedState() const
static constexpr int DT
Definition: MuonSubdetId.h:11
MuonTransientTrackingRecHit::MuonRecHitContainer unusedHits(const DetLayer *, const TrajectoryMeasurement &) const
double t0() const
Get the segment t0 (if recomputed, 0 is returned otherwise)
void reverseDirection(TrajectoryStateOnSurface &, const MagneticField *) const
void build(const TrajectoryStateOnSurface &, const NavigationDirection &, Trajectory &)
tmp
align.sh
Definition: createJobs.py:716
static constexpr int CSC
Definition: MuonSubdetId.h:12
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:53
std::vector< MuonRecHitPointer > MuonRecHitContainer
TrajectoryContainer trajectories(const Trajectory &traj) const override
Definition: event.py:1
GlobalVector globalDirection() const
T dot(const Basic3DVector &rh) const
Scalar product, or "dot" product, with a vector of same type.