CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GlobalCosmicMuonTrajectoryBuilder.cc
Go to the documentation of this file.
1 
14 
17 
26 
27 using namespace std;
28 using namespace edm;
29 
30 //
31 // constructor
32 //
34  const MuonServiceProxy* service) : theService(service) {
35  ParameterSet smootherPSet = par.getParameter<ParameterSet>("SmootherParameters");
36  theSmoother = new CosmicMuonSmoother(smootherPSet,theService);
37 
38  ParameterSet trackMatcherPSet = par.getParameter<ParameterSet>("GlobalMuonTrackMatcher");
39  theTrackMatcher = new GlobalMuonTrackMatcher(trackMatcherPSet,theService);
40 
41  theTkTrackLabel = par.getParameter<InputTag>("TkTrackCollectionLabel");
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 //
50 // destructor
51 //
52 
54 
55  if (theSmoother) delete theSmoother;
56  if (theTrackMatcher) delete theTrackMatcher;
57 }
58 
59 //
60 // set Event
61 //
63  event.getByLabel(theTkTrackLabel,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 //
83 // reconstruct trajectories
84 //
87 
88  if (!theTrackerTracks.isValid()) {
89  LogTrace(category_)<< "Tracker Track collection is invalid!!!";
90  return result;
91  }
92 
93  LogTrace(category_) <<"Found "<<theTrackerTracks->size()<<" tracker Tracks";
94  if (theTrackerTracks->empty()) 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()) return result;
101  reco::TrackRef tkTrack = matched.front().second;
102 
103  if ( tkTrack.isNull() ) return result;
104  reco::TrackRef muTrack = muCand.second;
105 
106  ConstRecHitContainer muRecHits;
107 
108  if (muCand.first == 0 || !muCand.first->isValid()) {
109  muRecHits = getTransientRecHits(*muTrack);
110  } else {
111  muRecHits = muCand.first->recHits();
112  }
113 
114  LogTrace(category_)<<"mu RecHits: "<<muRecHits.size();
115 
116  ConstRecHitContainer tkRecHits = getTransientRecHits(*tkTrack);
117 
118 // if ( !tkTrajsAvailable ) {
119 // tkRecHits = getTransientRecHits(*tkTrack);
120 // } else {
121 // tkRecHits = allTrackerTrajs->front().recHits();
122 // }
123 
124  ConstRecHitContainer hits; //= tkRecHits;
125  LogTrace(category_)<<"tk RecHits: "<<tkRecHits.size();
126 
127 // hits.insert(hits.end(), muRecHits.begin(), muRecHits.end());
128 // stable_sort(hits.begin(), hits.end(), DecreasingGlobalY());
129 
130  sortHits(hits, muRecHits, tkRecHits);
131 
132 // LogTrace(category_)<< "Used RecHits after sort: "<<hits.size()<<endl;;
133 // LogTrace(category_) <<utilities()->print(hits)<<endl;
134 // LogTrace(category_) << "== End of Used RecHits == "<<endl;
135 
136  TrajectoryStateTransform tsTrans;
137 
138  TrajectoryStateOnSurface muonState1 = tsTrans.innerStateOnSurface(*muTrack, *theService->trackingGeometry(), &*theService->magneticField());
139  TrajectoryStateOnSurface tkState1 = tsTrans.innerStateOnSurface(*tkTrack, *theService->trackingGeometry(), &*theService->magneticField());
140 
141  TrajectoryStateOnSurface muonState2 = tsTrans.outerStateOnSurface(*muTrack, *theService->trackingGeometry(), &*theService->magneticField());
142  TrajectoryStateOnSurface tkState2 = tsTrans.outerStateOnSurface(*tkTrack, *theService->trackingGeometry(), &*theService->magneticField());
143 
144  TrajectoryStateOnSurface firstState1 =
145  ( muonState1.globalPosition().y() > tkState1.globalPosition().y() )? muonState1 : tkState1;
146  TrajectoryStateOnSurface firstState2 =
147  ( muonState2.globalPosition().y() > tkState2.globalPosition().y() )? muonState2 : tkState2;
148 
149  TrajectoryStateOnSurface firstState =
150  ( firstState1.globalPosition().y() > firstState2.globalPosition().y() )? firstState1 : firstState2;
151 
153  if(hits.size()>0){
154  front=hits.front()->globalPosition();
155  back=hits.back()->globalPosition();
156  if ( (front.perp()<130 && fabs(front.z())<300)|| (back.perp() <130 && fabs(back.z())<300)){
157  if (hits.front()->globalPosition().perp()>hits.back()->globalPosition().perp()) reverse(hits.begin(), hits.end());
158  tkState1 = tsTrans.innerStateOnSurface(*tkTrack, *theService->trackingGeometry(), &*theService->magneticField());
159  tkState2 = tsTrans.outerStateOnSurface(*tkTrack, *theService->trackingGeometry(), &*theService->magneticField());
160  firstState =( tkState1.globalPosition().perp() < tkState2.globalPosition().perp() )? tkState1 : tkState2;
161  }
162  }
163  if (!firstState.isValid()) return result;
164  LogTrace(category_) <<"firstTSOS pos: "<<firstState.globalPosition()<<"mom: "<<firstState.globalMomentum();
165 
166  // begin refitting
167 
168  TrajectorySeed seed;
169  vector<Trajectory> refitted = theSmoother->trajectories(seed,hits,firstState);
170 
171  if ( refitted.empty() ) {
172  LogTrace(category_)<<"smoothing trajectories fail";
173  refitted = theSmoother->fit(seed,hits,firstState); //FIXME
174  }
175 
176  if (refitted.empty()) {
177  LogTrace(category_)<<"refit fail";
178  return result;
179  }
180 
181  Trajectory* myTraj = new Trajectory(refitted.front());
182 
183  const std::vector<TrajectoryMeasurement>& mytms = myTraj->measurements();
184  LogTrace(category_)<<"measurements in final trajectory "<<mytms.size();
185  LogTrace(category_) <<"Orignally there are "<<tkTrack->found()<<" tk rhs and "<<muTrack->found()<<" mu rhs.";
186 
187  if ( mytms.size() <= tkTrack->found() ) {
188  LogTrace(category_)<<"insufficient measurements. skip... ";
189  return result;
190  }
191 
192  MuonCandidate* myCand = new MuonCandidate(myTraj,muTrack,tkTrack);
193  result.push_back(myCand);
194  LogTrace(category_)<<"final global cosmic muon: ";
195  for (std::vector<TrajectoryMeasurement>::const_iterator itm = mytms.begin();
196  itm != mytms.end(); ++itm ) {
197  LogTrace(category_)<<"updated pos "<<itm->updatedState().globalPosition()
198  <<"mom "<<itm->updatedState().globalMomentum();
199  }
200  return result;
201 }
202 
204 
205  if ( tkHits.empty() ) {
206  LogTrace(category_) << "No valid tracker hits";
207  return;
208  }
209  if ( muonHits.empty() ) {
210  LogTrace(category_) << "No valid muon hits";
211  return;
212  }
213 
214  ConstRecHitContainer::const_iterator frontTkHit = tkHits.begin();
215  ConstRecHitContainer::const_iterator backTkHit = tkHits.end() - 1;
216  while ( !(*frontTkHit)->isValid() && frontTkHit != backTkHit) {frontTkHit++;}
217  while ( !(*backTkHit)->isValid() && backTkHit != frontTkHit) {backTkHit--;}
218 
219  ConstRecHitContainer::const_iterator frontMuHit = muonHits.begin();
220  ConstRecHitContainer::const_iterator backMuHit = muonHits.end() - 1;
221  while ( !(*frontMuHit)->isValid() && frontMuHit != backMuHit) {frontMuHit++;}
222  while ( !(*backMuHit)->isValid() && backMuHit != frontMuHit) {backMuHit--;}
223 
224  if ( frontTkHit == backTkHit ) {
225  LogTrace(category_) << "No valid tracker hits";
226  return;
227  }
228  if ( frontMuHit == backMuHit ) {
229  LogTrace(category_) << "No valid muon hits";
230  return;
231  }
232 
233  GlobalPoint frontTkPos = (*frontTkHit)->globalPosition();
234  GlobalPoint backTkPos = (*backTkHit)->globalPosition();
235 
236  GlobalPoint frontMuPos = (*frontMuHit)->globalPosition();
237  GlobalPoint backMuPos = (*backMuHit)->globalPosition();
238 
239  //sort hits going from higher to lower positions
240  if ( frontTkPos.y() < backTkPos.y() ) {//check if tk hits order same direction
241  reverse(tkHits.begin(), tkHits.end());
242  }
243 
244  if ( frontMuPos.y() < backMuPos.y() ) {
245  reverse(muonHits.begin(), muonHits.end());
246  }
247 
248 // LogTrace(category_)<< "tkHits after sort: "<<tkHits.size()<<endl;;
249 // LogTrace(category_) <<utilities()->print(tkHits)<<endl;
250 // LogTrace(category_) << "== End of tkHits == "<<endl;
251 
252 // LogTrace(category_)<< "muonHits after sort: "<<muonHits.size()<<endl;;
253 // LogTrace(category_) <<utilities()->print(muonHits)<<endl;
254 // LogTrace(category_)<< "== End of muonHits == "<<endl;
255 
256  //separate muon hits into 2 different hemisphere
257  ConstRecHitContainer::iterator middlepoint = muonHits.begin();
258  bool insertInMiddle = false;
259 
260  for (ConstRecHitContainer::iterator ihit = muonHits.begin();
261  ihit != muonHits.end() - 1; ihit++ ) {
262  GlobalPoint ipos = (*ihit)->globalPosition();
263  GlobalPoint nextpos = (*(ihit+1))->globalPosition();
264  if ( (ipos-nextpos).mag() < 100.0 ) continue;
265 
266  GlobalPoint middle((ipos.x()+nextpos.x())/2, (ipos.y()+nextpos.y())/2, (ipos.z()+nextpos.z())/2);
267  LogTrace(category_)<<"ipos "<<ipos<<"nextpos"<<nextpos<<" middle "<<middle<<endl;
268  if ( (middle.perp() < ipos.perp()) && (middle.perp() < nextpos.perp() ) ) {
269  LogTrace(category_)<<"found middlepoint"<<endl;
270  middlepoint = ihit;
271  insertInMiddle = true;
272  break;
273  }
274  }
275 
276  //insert track hits in correct order
277  if ( insertInMiddle ) { //if tk hits should be sandwich
278  GlobalPoint jointpointpos = (*middlepoint)->globalPosition();
279  LogTrace(category_)<<"jointpoint "<<jointpointpos<<endl;
280  if ((frontTkPos - jointpointpos).mag() > (backTkPos - jointpointpos).mag() ) {//check if tk hits order same direction
281  reverse(tkHits.begin(), tkHits.end());
282  }
283  muonHits.insert(middlepoint+1, tkHits.begin(), tkHits.end());
284  hits = muonHits;
285  } else { // append at one end
286  if ( frontTkPos.y() < frontMuPos.y() ) { //insert at the end
287  LogTrace(category_)<<"insert at the end "<<frontTkPos << frontMuPos <<endl;
288 
289  hits = muonHits;
290  hits.insert(hits.end(), tkHits.begin(), tkHits.end());
291  } else { //insert at the beginning
292  LogTrace(category_)<<"insert at the beginning "<<frontTkPos << frontMuPos <<endl;
293  hits = tkHits;
294  hits.insert(hits.end(), muonHits.begin(), muonHits.end());
295  }
296  }
297 }
298 
299 
302 
304 
305  TrajectoryStateTransform tsTrans;
306 
307  TrajectoryStateOnSurface currTsos = tsTrans.innerStateOnSurface(track, *theService->trackingGeometry(), &*theService->magneticField());
308  for (trackingRecHit_iterator hit = track.recHitsBegin(); hit != track.recHitsEnd(); ++hit) {
309  if((*hit)->isValid()) {
310  DetId recoid = (*hit)->geographicalId();
311  if ( recoid.det() == DetId::Tracker ) {
313  TrajectoryStateOnSurface predTsos = theService->propagator(thePropagatorName)->propagate(currTsos, theService->trackingGeometry()->idToDet(recoid)->surface());
314  LogTrace(category_)<<"predtsos "<<predTsos.isValid();
315  if ( predTsos.isValid() ) {
316  currTsos = predTsos;
317  TransientTrackingRecHit::RecHitPointer preciseHit = ttrhit->clone(currTsos);
318  result.push_back(preciseHit);
319  }
320  } else if ( recoid.det() == DetId::Muon ) {
321  result.push_back(theMuonRecHitBuilder->build(&**hit));
322  }
323  }
324  }
325  return result;
326 }
327 
328 std::vector<GlobalCosmicMuonTrajectoryBuilder::TrackCand> GlobalCosmicMuonTrajectoryBuilder::match(const TrackCand& mu, const edm::Handle<reco::TrackCollection>& tktracks) {
329 
330  std::vector<TrackCand> result;
331 
332  TrajectoryStateTransform tsTrans;
333  TrajectoryStateOnSurface innerTsos = tsTrans.innerStateOnSurface(*(mu.second), *theService->trackingGeometry(), &*theService->magneticField());
334  TrajectoryStateOnSurface outerTsos = tsTrans.outerStateOnSurface(*(mu.second), *theService->trackingGeometry(), &*theService->magneticField());
335  //build tracker TrackCands and pick the best match if size greater than 2
336  vector<TrackCand> tkTrackCands;
339  TrackCand tkCand = TrackCand((Trajectory*)(0),tkTrack);
340  tkTrackCands.push_back(tkCand);
341  LogTrace(category_) << "chisq is " << theTrackMatcher->match(mu,tkCand,0,0);
342  LogTrace(category_) << "d is " << theTrackMatcher->match(mu,tkCand,1,0);
343  LogTrace(category_) << "r_pos is " << theTrackMatcher->match(mu,tkCand,2,0);
344  }
345 
346  // now if only 1 tracker tracks, return it
347  if (tkTrackCands.size() <= 1 ) {
348  return tkTrackCands;
349  }
350 
351  // if there're many tracker tracks
352 
353  // if muon is only on one side
354  GlobalPoint innerPos = innerTsos.globalPosition();
355  GlobalPoint outerPos = outerTsos.globalPosition();
356 
357  if ( ( innerPos.basicVector().dot( innerTsos.globalMomentum().basicVector() ) *
358  outerPos.basicVector().dot(outerTsos.globalMomentum().basicVector() ) > 0 ) ) {
359 
360  GlobalPoint geoInnerPos = (innerPos.mag() < outerPos.mag()) ? innerPos : outerPos;
361  LogTrace(category_) <<"geoInnerPos Mu "<<geoInnerPos<<endl;
362 
363  // if there're tracker tracks totally on the other half
364  // and there're tracker tracks on the same half
365  // remove the tracks on the other half
366  for(vector<TrackCand>::const_iterator itkCand = tkTrackCands.begin(); itkCand != tkTrackCands.end(); ++itkCand) {
367 
368  reco::TrackRef tkTrack = itkCand->second;
369 
370  GlobalPoint tkInnerPos(tkTrack->innerPosition().x(), tkTrack->innerPosition().y(), tkTrack->innerPosition().z());
371  GlobalPoint tkOuterPos(tkTrack->outerPosition().x(), tkTrack->outerPosition().y(), tkTrack->outerPosition().z());
372  LogTrace(category_) <<"tkTrack "<<tkInnerPos<<" "<<tkOuterPos<<endl;
373 
374  float closetDistance11 = (geoInnerPos - tkInnerPos).mag() ;
375  float closetDistance12 = (geoInnerPos - tkOuterPos).mag() ;
376  float closetDistance1 = (closetDistance11 < closetDistance12) ? closetDistance11 : closetDistance12;
377  LogTrace(category_) <<"closetDistance1 "<<closetDistance1<<endl;
378 
379  if (true || !isTraversing(*tkTrack) ) {
380  bool keep = true;
381  for(vector<TrackCand>::const_iterator itkCand2 = tkTrackCands.begin(); itkCand2 != tkTrackCands.end(); ++itkCand2) {
382  if (itkCand2 == itkCand ) continue;
383  reco::TrackRef tkTrack2 = itkCand2->second;
384 
385  GlobalPoint tkInnerPos2(tkTrack2->innerPosition().x(), tkTrack2->innerPosition().y(), tkTrack2->innerPosition().z());
386  GlobalPoint tkOuterPos2(tkTrack2->outerPosition().x(), tkTrack2->outerPosition().y(), tkTrack2->outerPosition().z());
387  LogTrace(category_) <<"tkTrack2 "<< tkInnerPos2 <<" "<<tkOuterPos2 <<endl;
388 
389  float farthestDistance21 = (geoInnerPos - tkInnerPos2).mag() ;
390  float farthestDistance22 = (geoInnerPos - tkOuterPos2).mag() ;
391  float farthestDistance2 = (farthestDistance21 > farthestDistance22) ? farthestDistance21 : farthestDistance22;
392  LogTrace(category_) <<"farthestDistance2 "<<farthestDistance2<<endl;
393 
394  if (closetDistance1 > farthestDistance2 - 1e-3) {
395  keep = false;
396  break;
397  }
398  }
399  if (keep) result.push_back(*itkCand);
400  else LogTrace(category_) <<"The Track is on different hemisphere"<<endl;
401  } else {
402  result.push_back(*itkCand);
403  }
404  }
405  if ( result.empty() ) {
406  //if all tk tracks on the other side, still keep them
407  result = tkTrackCands;
408  }
409  } else { // muon is traversing
410  result = tkTrackCands;
411  }
412 
413  // match muCand to tkTrackCands
414  vector<TrackCand> matched_trackerTracks = theTrackMatcher->match(mu, result);
415 
416  LogTrace(category_) <<"TrackMatcher found " << matched_trackerTracks.size() << "tracker tracks matched";
417 
418  //now pick the best matched one
419  if( matched_trackerTracks.size() < 2 ) {
420  return matched_trackerTracks;
421  } else {
422  // in case of more than 1 tkTrack,
423  // select the best-one based on distance (matchOption==1)
424  // at innermost Mu hit surface. (surfaceOption == 0)
425  result.clear();
426  TrackCand bestMatch;
427 
428  double quality = 1e6;
429  double max_quality = 1e6;
430  for( vector<TrackCand>::const_iterator iter = matched_trackerTracks.begin(); iter != matched_trackerTracks.end(); iter++) {
431  quality = theTrackMatcher->match(mu,*iter, 1, 0);
432  LogTrace(category_) <<" quality of tracker track is " << quality;
433  if( quality < max_quality ) {
434  max_quality=quality;
435  bestMatch = (*iter);
436  }
437  }
438  LogTrace(category_) <<" Picked tracker track with quality " << max_quality;
439  result.push_back(bestMatch);
440  return result;
441  }
442 }
443 
444 
446 
448  for (trackingRecHit_iterator hit = track.recHitsBegin(); hit != track.recHitsEnd(); ++hit) {
449  if((*hit)->isValid()) {
450  firstValid = hit;
451  break;
452  }
453  }
454 
455  trackingRecHit_iterator lastValid;
456  for (trackingRecHit_iterator hit = track.recHitsEnd() - 1; hit != track.recHitsBegin() - 1; --hit) {
457  if((*hit)->isValid()) {
458  lastValid = hit;
459  break;
460  }
461  }
462 
463  GlobalPoint posFirst = theService->trackingGeometry()->idToDet((*firstValid)->geographicalId())->position();
464 
465  GlobalPoint posLast = theService->trackingGeometry()->idToDet((*lastValid)->geographicalId() )->position();
466 
467  GlobalPoint middle((posFirst.x()+posLast.x())/2, (posFirst.y()+posLast.y())/2, (posFirst.z()+posLast.z())/2);
468 
469  if ( (middle.mag() < posFirst.mag()) && (middle.mag() < posLast.mag() ) ) {
470  return true;
471  }
472  return false;
473 
474 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
T perp() const
Definition: PV3DBase.h:66
edm::ESHandle< TransientTrackingRecHitBuilder > theTrackerRecHitBuilder
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
edm::ESHandle< TransientTrackingRecHitBuilder > theMuonRecHitBuilder
T y() const
Definition: PV3DBase.h:57
GlobalPoint globalPosition() const
std::vector< Trajectory * > trajectories(const TrajectorySeed &)
dummy implementation, unused in this class
uint16_t size_type
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
TransientTrackingRecHit::ConstRecHitContainer ConstRecHitContainer
DataContainer const & measurements() const
Definition: Trajectory.h:169
const int keep
ConstRecHitContainer getTransientRecHits(const reco::Track &) const
bool isNull() const
Checks for null.
Definition: Ref.h:246
T mag() const
Definition: PV3DBase.h:61
std::pair< const Trajectory *, reco::TrackRef > TrackCand
bool isTraversing(const reco::Track &tk) const
T z() const
Definition: PV3DBase.h:58
tuple result
Definition: query.py:137
TrajectoryStateOnSurface outerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field) const
void sortHits(ConstRecHitContainer &, ConstRecHitContainer &, ConstRecHitContainer &)
virtual void setEvent(const edm::Event &)
pass the Event to the algo at each event
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:63
std::vector< MuonCandidate * > CandidateContainer
Definition: MuonCandidate.h:22
double match(const TrackCand &sta, const TrackCand &track, int matchOption=0, int surfaceOption=1) const
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool isValid() const
Definition: HandleBase.h:76
#define LogTrace(id)
std::vector< ConstRecHitPointer > ConstRecHitContainer
Definition: DetId.h:20
std::vector< Trajectory > fit(const Trajectory &) const
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field) const
std::vector< TrackCand > match(const TrackCand &, const edm::Handle< reco::TrackCollection > &)
check if tk and muon Tracks are matched
GlobalVector globalMomentum() const
virtual std::vector< Trajectory > trajectories(const Trajectory &) const
Detector det() const
get the detector field from this detid
Definition: DetId.h:37
T x() const
Definition: PV3DBase.h:56
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:54
edm::Handle< reco::TrackCollection > theTrackerTracks
GlobalCosmicMuonTrajectoryBuilder(const edm::ParameterSet &, const MuonServiceProxy *service)
Constructor.
const double par[8 *NPar][4]
T dot(const Basic3DVector &rh) const
Scalar product, or &quot;dot&quot; product, with a vector of same type.
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:65