CMS 3D CMS Logo

GlobalCosmicMuonTrajectoryBuilder.cc
Go to the documentation of this file.
1 
11 
14 
23 
24 using namespace std;
25 using namespace edm;
26 
27 //
28 // constructor
29 //
33  : theService(service) {
34  ParameterSet smootherPSet = par.getParameter<ParameterSet>("SmootherParameters");
35  theSmoother = new CosmicMuonSmoother(smootherPSet, theService);
36 
37  ParameterSet trackMatcherPSet = par.getParameter<ParameterSet>("GlobalMuonTrackMatcher");
38  theTrackMatcher = new GlobalMuonTrackMatcher(trackMatcherPSet, theService);
39 
40  theTkTrackToken = iC.consumes<reco::TrackCollection>(par.getParameter<InputTag>("TkTrackCollectionLabel"));
41 
42  theTrackerRecHitBuilderName = par.getParameter<string>("TrackerRecHitBuilder");
43  theMuonRecHitBuilderName = par.getParameter<string>("MuonRecHitBuilder");
44  thePropagatorName = par.getParameter<string>("Propagator");
45  category_ = "Muon|RecoMuon|CosmicMuon|GlobalCosmicMuonTrajectoryBuilder";
46 }
47 
48 //
49 // destructor
50 //
51 
53  if (theSmoother)
54  delete theSmoother;
55  if (theTrackMatcher)
56  delete theTrackMatcher;
57 }
58 
59 //
60 // set Event
61 //
63  event.getByToken(theTkTrackToken, theTrackerTracks);
64 
65  // edm::Handle<std::vector<Trajectory> > handleTrackerTrajs;
66  // if ( event.getByLabel(theTkTrackLabel,handleTrackerTrajs) && handleTrackerTrajs.isValid() ) {
67  // tkTrajsAvailable = true;
68  // allTrackerTrajs = &*handleTrackerTrajs;
69  // LogInfo("GlobalCosmicMuonTrajectoryBuilder")
70  // << "Tk Trajectories Found! " << endl;
71  // } else {
72  // LogInfo("GlobalCosmicMuonTrajectoryBuilder")
73  // << "No Tk Trajectories Found! " << endl;
74  // tkTrajsAvailable = false;
75  // }
76 
79 }
80 
81 //
82 // reconstruct trajectories
83 //
86 
87  if (!theTrackerTracks.isValid()) {
88  LogTrace(category_) << "Tracker Track collection is invalid!!!";
89  return result;
90  }
91 
92  LogTrace(category_) << "Found " << theTrackerTracks->size() << " tracker Tracks";
93  if (theTrackerTracks->empty())
94  return result;
95 
96  vector<TrackCand> matched = match(muCand, theTrackerTracks);
97 
98  LogTrace(category_) << "TrackMatcher found " << matched.size() << "tracker tracks matched";
99 
100  if (matched.empty())
101  return result;
102  reco::TrackRef tkTrack = matched.front().second;
103 
104  if (tkTrack.isNull())
105  return result;
106  reco::TrackRef muTrack = muCand.second;
107 
108  ConstRecHitContainer muRecHits;
109 
110  if (muCand.first == nullptr || !muCand.first->isValid()) {
111  muRecHits = getTransientRecHits(*muTrack);
112  } else {
113  muRecHits = muCand.first->recHits();
114  }
115 
116  LogTrace(category_) << "mu RecHits: " << muRecHits.size();
117 
118  ConstRecHitContainer tkRecHits = getTransientRecHits(*tkTrack);
119 
120  // if ( !tkTrajsAvailable ) {
121  // tkRecHits = getTransientRecHits(*tkTrack);
122  // } else {
123  // tkRecHits = allTrackerTrajs->front().recHits();
124  // }
125 
126  ConstRecHitContainer hits; //= tkRecHits;
127  LogTrace(category_) << "tk RecHits: " << tkRecHits.size();
128 
129  // hits.insert(hits.end(), muRecHits.begin(), muRecHits.end());
130  // stable_sort(hits.begin(), hits.end(), DecreasingGlobalY());
131 
132  sortHits(hits, muRecHits, tkRecHits);
133 
134  // LogTrace(category_)<< "Used RecHits after sort: "<<hits.size()<<endl;;
135  // LogTrace(category_) <<utilities()->print(hits)<<endl;
136  // LogTrace(category_) << "== End of Used RecHits == "<<endl;
137 
139  *muTrack, *theService->trackingGeometry(), &*theService->magneticField());
141  *tkTrack, *theService->trackingGeometry(), &*theService->magneticField());
142 
144  *muTrack, *theService->trackingGeometry(), &*theService->magneticField());
146  *tkTrack, *theService->trackingGeometry(), &*theService->magneticField());
147 
148  TrajectoryStateOnSurface firstState1 =
149  (muonState1.globalPosition().y() > tkState1.globalPosition().y()) ? muonState1 : tkState1;
150  TrajectoryStateOnSurface firstState2 =
151  (muonState2.globalPosition().y() > tkState2.globalPosition().y()) ? muonState2 : tkState2;
152 
153  TrajectoryStateOnSurface firstState =
154  (firstState1.globalPosition().y() > firstState2.globalPosition().y()) ? firstState1 : firstState2;
155 
156  GlobalPoint front, back;
157  if (!hits.empty()) {
158  front = hits.front()->globalPosition();
159  back = hits.back()->globalPosition();
160  if ((front.perp() < 130 && fabs(front.z()) < 300) || (back.perp() < 130 && fabs(back.z()) < 300)) {
161  if (hits.front()->globalPosition().perp() > hits.back()->globalPosition().perp())
162  reverse(hits.begin(), hits.end());
164  *tkTrack, *theService->trackingGeometry(), &*theService->magneticField());
166  *tkTrack, *theService->trackingGeometry(), &*theService->magneticField());
167  firstState = (tkState1.globalPosition().perp() < tkState2.globalPosition().perp()) ? tkState1 : tkState2;
168  }
169  }
170  if (!firstState.isValid())
171  return result;
172  LogTrace(category_) << "firstTSOS pos: " << firstState.globalPosition() << "mom: " << firstState.globalMomentum();
173 
174  // begin refitting
175 
177  vector<Trajectory> refitted = theSmoother->trajectories(seed, hits, firstState);
178 
179  if (refitted.empty()) {
180  LogTrace(category_) << "smoothing trajectories fail";
181  refitted = theSmoother->fit(seed, hits, firstState); //FIXME
182  }
183 
184  if (refitted.empty()) {
185  LogTrace(category_) << "refit fail";
186  return result;
187  }
188 
189  Trajectory* myTraj = new Trajectory(refitted.front());
190 
191  const std::vector<TrajectoryMeasurement>& mytms = myTraj->measurements();
192  LogTrace(category_) << "measurements in final trajectory " << mytms.size();
193  LogTrace(category_) << "Orignally there are " << tkTrack->found() << " tk rhs and " << muTrack->found() << " mu rhs.";
194 
195  if (mytms.size() <= tkTrack->found()) {
196  LogTrace(category_) << "insufficient measurements. skip... ";
197  return result;
198  }
199 
200  MuonCandidate* myCand = new MuonCandidate(myTraj, muTrack, tkTrack);
201  result.push_back(myCand);
202  LogTrace(category_) << "final global cosmic muon: ";
203  for (std::vector<TrajectoryMeasurement>::const_iterator itm = mytms.begin(); itm != mytms.end(); ++itm) {
204  LogTrace(category_) << "updated pos " << itm->updatedState().globalPosition() << "mom "
205  << itm->updatedState().globalMomentum();
206  }
207  return result;
208 }
209 
211  ConstRecHitContainer& muonHits,
212  ConstRecHitContainer& tkHits) {
213  if (tkHits.empty()) {
214  LogTrace(category_) << "No valid tracker hits";
215  return;
216  }
217  if (muonHits.empty()) {
218  LogTrace(category_) << "No valid muon hits";
219  return;
220  }
221 
222  ConstRecHitContainer::const_iterator frontTkHit = tkHits.begin();
223  ConstRecHitContainer::const_iterator backTkHit = tkHits.end() - 1;
224  while (!(*frontTkHit)->isValid() && frontTkHit != backTkHit) {
225  frontTkHit++;
226  }
227  while (!(*backTkHit)->isValid() && backTkHit != frontTkHit) {
228  backTkHit--;
229  }
230 
231  ConstRecHitContainer::const_iterator frontMuHit = muonHits.begin();
232  ConstRecHitContainer::const_iterator backMuHit = muonHits.end() - 1;
233  while (!(*frontMuHit)->isValid() && frontMuHit != backMuHit) {
234  frontMuHit++;
235  }
236  while (!(*backMuHit)->isValid() && backMuHit != frontMuHit) {
237  backMuHit--;
238  }
239 
240  if (frontTkHit == backTkHit) {
241  LogTrace(category_) << "No valid tracker hits";
242  return;
243  }
244  if (frontMuHit == backMuHit) {
245  LogTrace(category_) << "No valid muon hits";
246  return;
247  }
248 
249  GlobalPoint frontTkPos = (*frontTkHit)->globalPosition();
250  GlobalPoint backTkPos = (*backTkHit)->globalPosition();
251 
252  GlobalPoint frontMuPos = (*frontMuHit)->globalPosition();
253  GlobalPoint backMuPos = (*backMuHit)->globalPosition();
254 
255  //sort hits going from higher to lower positions
256  if (frontTkPos.y() < backTkPos.y()) { //check if tk hits order same direction
257  reverse(tkHits.begin(), tkHits.end());
258  }
259 
260  if (frontMuPos.y() < backMuPos.y()) {
261  reverse(muonHits.begin(), muonHits.end());
262  }
263 
264  // LogTrace(category_)<< "tkHits after sort: "<<tkHits.size()<<endl;;
265  // LogTrace(category_) <<utilities()->print(tkHits)<<endl;
266  // LogTrace(category_) << "== End of tkHits == "<<endl;
267 
268  // LogTrace(category_)<< "muonHits after sort: "<<muonHits.size()<<endl;;
269  // LogTrace(category_) <<utilities()->print(muonHits)<<endl;
270  // LogTrace(category_)<< "== End of muonHits == "<<endl;
271 
272  //separate muon hits into 2 different hemisphere
273  ConstRecHitContainer::iterator middlepoint = muonHits.begin();
274  bool insertInMiddle = false;
275 
276  for (ConstRecHitContainer::iterator ihit = muonHits.begin(); ihit != muonHits.end() - 1; ihit++) {
277  GlobalPoint ipos = (*ihit)->globalPosition();
278  GlobalPoint nextpos = (*(ihit + 1))->globalPosition();
279  if ((ipos - nextpos).mag() < 100.0)
280  continue;
281 
282  GlobalPoint middle((ipos.x() + nextpos.x()) / 2, (ipos.y() + nextpos.y()) / 2, (ipos.z() + nextpos.z()) / 2);
283  LogTrace(category_) << "ipos " << ipos << "nextpos" << nextpos << " middle " << middle << endl;
284  if ((middle.perp() < ipos.perp()) && (middle.perp() < nextpos.perp())) {
285  LogTrace(category_) << "found middlepoint" << endl;
286  middlepoint = ihit;
287  insertInMiddle = true;
288  break;
289  }
290  }
291 
292  //insert track hits in correct order
293  if (insertInMiddle) { //if tk hits should be sandwich
294  GlobalPoint jointpointpos = (*middlepoint)->globalPosition();
295  LogTrace(category_) << "jointpoint " << jointpointpos << endl;
296  if ((frontTkPos - jointpointpos).mag() >
297  (backTkPos - jointpointpos).mag()) { //check if tk hits order same direction
298  reverse(tkHits.begin(), tkHits.end());
299  }
300  muonHits.insert(middlepoint + 1, tkHits.begin(), tkHits.end());
301  hits = muonHits;
302  } else { // append at one end
303  if (frontTkPos.y() < frontMuPos.y()) { //insert at the end
304  LogTrace(category_) << "insert at the end " << frontTkPos << frontMuPos << endl;
305 
306  hits = muonHits;
307  hits.insert(hits.end(), tkHits.begin(), tkHits.end());
308  } else { //insert at the beginning
309  LogTrace(category_) << "insert at the beginning " << frontTkPos << frontMuPos << endl;
310  hits = tkHits;
311  hits.insert(hits.end(), muonHits.begin(), muonHits.end());
312  }
313  }
314 }
315 
317 
319  const reco::Track& track) const {
321 
322  auto hitCloner = static_cast<TkTransientTrackingRecHitBuilder const*>(theTrackerRecHitBuilder.product())->cloner();
324  track, *theService->trackingGeometry(), &*theService->magneticField());
325  for (trackingRecHit_iterator hit = track.recHitsBegin(); hit != track.recHitsEnd(); ++hit) {
326  if ((*hit)->isValid()) {
327  DetId recoid = (*hit)->geographicalId();
328  if (recoid.det() == DetId::Tracker) {
329  TrajectoryStateOnSurface predTsos =
330  theService->propagator(thePropagatorName)
331  ->propagate(currTsos, theService->trackingGeometry()->idToDet(recoid)->surface());
332  LogTrace(category_) << "predtsos " << predTsos.isValid();
333  if (predTsos.isValid()) {
334  currTsos = predTsos;
335  result.emplace_back(hitCloner(**hit, predTsos));
336  }
337  } else if (recoid.det() == DetId::Muon) {
338  result.push_back(theMuonRecHitBuilder->build(&**hit));
339  }
340  }
341  }
342  return result;
343 }
344 
345 std::vector<GlobalCosmicMuonTrajectoryBuilder::TrackCand> GlobalCosmicMuonTrajectoryBuilder::match(
346  const TrackCand& mu, const edm::Handle<reco::TrackCollection>& tktracks) {
347  std::vector<TrackCand> result;
348 
350  *(mu.second), *theService->trackingGeometry(), &*theService->magneticField());
352  *(mu.second), *theService->trackingGeometry(), &*theService->magneticField());
353  //build tracker TrackCands and pick the best match if size greater than 2
354  vector<TrackCand> tkTrackCands;
355  for (reco::TrackCollection::size_type i = 0; i < theTrackerTracks->size(); ++i) {
357  TrackCand tkCand = TrackCand((Trajectory*)nullptr, tkTrack);
358  tkTrackCands.push_back(tkCand);
359  LogTrace(category_) << "chisq is " << theTrackMatcher->match(mu, tkCand, 0, 0);
360  LogTrace(category_) << "d is " << theTrackMatcher->match(mu, tkCand, 1, 0);
361  LogTrace(category_) << "r_pos is " << theTrackMatcher->match(mu, tkCand, 2, 0);
362  }
363 
364  // now if only 1 tracker tracks, return it
365  if (tkTrackCands.size() <= 1) {
366  return tkTrackCands;
367  }
368 
369  // if there're many tracker tracks
370 
371  // if muon is only on one side
372  GlobalPoint innerPos = innerTsos.globalPosition();
373  GlobalPoint outerPos = outerTsos.globalPosition();
374 
375  if ((innerPos.basicVector().dot(innerTsos.globalMomentum().basicVector()) *
376  outerPos.basicVector().dot(outerTsos.globalMomentum().basicVector()) >
377  0)) {
378  GlobalPoint geoInnerPos = (innerPos.mag() < outerPos.mag()) ? innerPos : outerPos;
379  LogTrace(category_) << "geoInnerPos Mu " << geoInnerPos << endl;
380 
381  // if there're tracker tracks totally on the other half
382  // and there're tracker tracks on the same half
383  // remove the tracks on the other half
384  for (vector<TrackCand>::const_iterator itkCand = tkTrackCands.begin(); itkCand != tkTrackCands.end(); ++itkCand) {
385  reco::TrackRef tkTrack = itkCand->second;
386 
387  GlobalPoint tkInnerPos(tkTrack->innerPosition().x(), tkTrack->innerPosition().y(), tkTrack->innerPosition().z());
388  GlobalPoint tkOuterPos(tkTrack->outerPosition().x(), tkTrack->outerPosition().y(), tkTrack->outerPosition().z());
389  LogTrace(category_) << "tkTrack " << tkInnerPos << " " << tkOuterPos << endl;
390 
391  float closetDistance11 = (geoInnerPos - tkInnerPos).mag();
392  float closetDistance12 = (geoInnerPos - tkOuterPos).mag();
393  float closetDistance1 = (closetDistance11 < closetDistance12) ? closetDistance11 : closetDistance12;
394  LogTrace(category_) << "closetDistance1 " << closetDistance1 << endl;
395 
396  if (true || !isTraversing(*tkTrack)) {
397  bool keep = true;
398  for (vector<TrackCand>::const_iterator itkCand2 = tkTrackCands.begin(); itkCand2 != tkTrackCands.end();
399  ++itkCand2) {
400  if (itkCand2 == itkCand)
401  continue;
402  reco::TrackRef tkTrack2 = itkCand2->second;
403 
404  GlobalPoint tkInnerPos2(
405  tkTrack2->innerPosition().x(), tkTrack2->innerPosition().y(), tkTrack2->innerPosition().z());
406  GlobalPoint tkOuterPos2(
407  tkTrack2->outerPosition().x(), tkTrack2->outerPosition().y(), tkTrack2->outerPosition().z());
408  LogTrace(category_) << "tkTrack2 " << tkInnerPos2 << " " << tkOuterPos2 << endl;
409 
410  float farthestDistance21 = (geoInnerPos - tkInnerPos2).mag();
411  float farthestDistance22 = (geoInnerPos - tkOuterPos2).mag();
412  float farthestDistance2 = (farthestDistance21 > farthestDistance22) ? farthestDistance21 : farthestDistance22;
413  LogTrace(category_) << "farthestDistance2 " << farthestDistance2 << endl;
414 
415  if (closetDistance1 > farthestDistance2 - 1e-3) {
416  keep = false;
417  break;
418  }
419  }
420  if (keep)
421  result.push_back(*itkCand);
422  else
423  LogTrace(category_) << "The Track is on different hemisphere" << endl;
424  } else {
425  result.push_back(*itkCand);
426  }
427  }
428  if (result.empty()) {
429  //if all tk tracks on the other side, still keep them
430  result = tkTrackCands;
431  }
432  } else { // muon is traversing
433  result = tkTrackCands;
434  }
435 
436  // match muCand to tkTrackCands
437  vector<TrackCand> matched_trackerTracks = theTrackMatcher->match(mu, result);
438 
439  LogTrace(category_) << "TrackMatcher found " << matched_trackerTracks.size() << "tracker tracks matched";
440 
441  //now pick the best matched one
442  if (matched_trackerTracks.size() < 2) {
443  return matched_trackerTracks;
444  } else {
445  // in case of more than 1 tkTrack,
446  // select the best-one based on distance (matchOption==1)
447  // at innermost Mu hit surface. (surfaceOption == 0)
448  result.clear();
450 
451  double quality = 1e6;
452  double max_quality = 1e6;
453  for (vector<TrackCand>::const_iterator iter = matched_trackerTracks.begin(); iter != matched_trackerTracks.end();
454  iter++) {
455  quality = theTrackMatcher->match(mu, *iter, 1, 0);
456  LogTrace(category_) << " quality of tracker track is " << quality;
457  if (quality < max_quality) {
458  max_quality = quality;
459  bestMatch = (*iter);
460  }
461  }
462  LogTrace(category_) << " Picked tracker track with quality " << max_quality;
463  result.push_back(bestMatch);
464  return result;
465  }
466 }
467 
470  for (trackingRecHit_iterator hit = track.recHitsBegin(); hit != track.recHitsEnd(); ++hit) {
471  if ((*hit)->isValid()) {
472  firstValid = hit;
473  break;
474  }
475  }
476 
477  trackingRecHit_iterator lastValid;
478  for (trackingRecHit_iterator hit = track.recHitsEnd() - 1; hit != track.recHitsBegin() - 1; --hit) {
479  if ((*hit)->isValid()) {
480  lastValid = hit;
481  break;
482  }
483  }
484 
485  GlobalPoint posFirst = theService->trackingGeometry()->idToDet((*firstValid)->geographicalId())->position();
486 
487  GlobalPoint posLast = theService->trackingGeometry()->idToDet((*lastValid)->geographicalId())->position();
488 
489  GlobalPoint middle(
490  (posFirst.x() + posLast.x()) / 2, (posFirst.y() + posLast.y()) / 2, (posFirst.z() + posLast.z()) / 2);
491 
492  if ((middle.mag() < posFirst.mag()) && (middle.mag() < posLast.mag())) {
493  return true;
494  }
495  return false;
496 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
T perp() const
Definition: PV3DBase.h:69
edm::ESHandle< TransientTrackingRecHitBuilder > theTrackerRecHitBuilder
edm::EDGetTokenT< reco::TrackCollection > theTkTrackToken
TrajectoryStateOnSurface outerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field, bool withErr=true)
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
edm::ESHandle< TransientTrackingRecHitBuilder > theMuonRecHitBuilder
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
T y() const
Definition: PV3DBase.h:60
GlobalPoint globalPosition() const
uint16_t size_type
GlobalCosmicMuonTrajectoryBuilder(const edm::ParameterSet &, const MuonServiceProxy *service, edm::ConsumesCollector &iC)
Constructor.
TransientTrackingRecHit::ConstRecHitContainer ConstRecHitContainer
DataContainer const & measurements() const
Definition: Trajectory.h:178
const int keep
ConstRecHitContainer getTransientRecHits(const reco::Track &) const
T mag() const
Definition: PV3DBase.h:64
bool isTraversing(const reco::Track &tk) const
std::vector< Trajectory * > trajectories(const TrajectorySeed &) override
dummy implementation, unused in this class
T z() const
Definition: PV3DBase.h:61
void sortHits(ConstRecHitContainer &, ConstRecHitContainer &, ConstRecHitContainer &)
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:88
std::vector< MuonCandidate * > CandidateContainer
Definition: MuonCandidate.h:17
double match(const TrackCand &sta, const TrackCand &track, int matchOption=0, int surfaceOption=1) const
bool isValid() const
Definition: HandleBase.h:70
bool isNull() const
Checks for null.
Definition: Ref.h:235
#define LogTrace(id)
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
std::vector< ConstRecHitPointer > ConstRecHitContainer
Definition: DetId.h:17
std::vector< Trajectory > fit(const Trajectory &) const
def bestMatch(object, matchCollection)
Definition: deltar.py:138
std::vector< TrackCand > match(const TrackCand &, const edm::Handle< reco::TrackCollection > &)
check if tk and muon Tracks are matched
void setEvent(const edm::Event &) override
pass the Event to the algo at each event
HLT enums.
GlobalVector globalMomentum() const
static int position[264][3]
Definition: ReadPGInfo.cc:289
std::pair< const Trajectory *, reco::TrackRef > TrackCand
T x() const
Definition: PV3DBase.h:59
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:53
T const * product() const
Definition: ESHandle.h:86
TrajectoryContainer trajectories(const Trajectory &traj) const override
edm::Handle< reco::TrackCollection > theTrackerTracks
Definition: event.py:1
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field, bool withErr=true)
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
T dot(const Basic3DVector &rh) const
Scalar product, or "dot" product, with a vector of same type.
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:91