CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackMerger.cc
Go to the documentation of this file.
1 #include "TrackMerger.h"
4 
5 
7 
10 
11 #include <boost/foreach.hpp>
12 #define foreach BOOST_FOREACH
13 
14 #define TRACK_SORT 1 // just use all hits from inner track, then append from outer outside it
15 #define DET_SORT 0 // sort hits using global position of detector
16 #define HIT_SORT 0 // sort hits using global position of hit
17 
18 
20  useInnermostState_(iConfig.getParameter<bool>("useInnermostState")),
21  debug_(iConfig.getUntrackedParameter<bool>("debug",false)),
22  theBuilderName(iConfig.getParameter<std::string>("ttrhBuilderName"))
23 {
24 }
25 
27 {
28 }
29 
30 void TrackMerger::init(const edm::EventSetup &iSetup)
31 {
35  iSetup.get<IdealGeometryRecord>().get(theTrkTopo);
36 }
37 
39 {
40  std::vector<const TrackingRecHit *> hits;
41  hits.reserve(inner.recHitsSize() + outer.recHitsSize());
42  if (debug_) std::cout << "Inner track hits: " << std::endl;
43  for (trackingRecHit_iterator it = inner.recHitsBegin(), ed = inner.recHitsEnd(); it != ed; ++it) {
44  hits.push_back(&**it);
45  if (debug_) {
46  DetId id(hits.back()->geographicalId());
47  std::cout << " subdet " << id.subdetId() << " layer " << theTrkTopo->layer(id) << " valid " << hits.back()->isValid() << " detid: " << id() << std::endl;
48  }
49  }
50  if (debug_) std::cout << "Outer track hits: " << std::endl;
51 
52 #if TRACK_SORT
53  DetId lastId(hits.back()->geographicalId());
54  int lastSubdet = lastId.subdetId();
55  unsigned int lastLayer = theTrkTopo->layer(lastId);
56  for (trackingRecHit_iterator it = outer.recHitsBegin(), ed = outer.recHitsEnd(); it != ed; ++it) {
57  const TrackingRecHit *hit = &**it;
58  DetId id(hit->geographicalId());
59  int thisSubdet = id.subdetId();
60  if (thisSubdet > lastSubdet || (thisSubdet == lastSubdet && theTrkTopo->layer(id) > lastLayer)) {
61  hits.push_back(hit);
62  if (debug_) std::cout << " adding subdet " << id.subdetId() << " layer " << theTrkTopo->layer(id) << " valid " << hit->isValid() << " detid: " << id() << std::endl;
63  } else {
64  if (debug_) std::cout << " skipping subdet " << thisSubdet << " layer " << theTrkTopo->layer(id) << " valid " << hit->isValid() << " detid: " << id() << std::endl;
65  }
66  }
67 #else
68  size_t nHitsFirstTrack = hits.size();
69  for (trackingRecHit_iterator it = outer.recHitsBegin(), ed = outer.recHitsEnd(); it != ed; ++it) {
70  const TrackingRecHit *hit = &**it;
71  DetId id(hit->geographicalId());
72  int lay = theTrkTopo->layer(id);
73  bool shared = false;
74  bool valid = hit->isValid();
75  if (debug_) std::cout << " subdet " << id.subdetId() << " layer " << theTrkTopo->layer(id) << " valid " << valid << " detid: " << id() << std::endl;
76  size_t iHit = 0;
77  foreach(const TrackingRecHit *& hit2, hits) {
78  ++iHit; if (iHit > nHitsFirstTrack) break;
79  DetId id2 = hit2->geographicalId();
80  if (id.subdetId() != id2.subdetId()) continue;
81  if (theTrkTopo->layer(id2) != lay) continue;
82  if (hit->sharesInput(hit2, TrackingRecHit::all)) {
83  if (debug_) std::cout << " discared as duplicate of other hit" << id() << std::endl;
84  shared = true; break;
85  }
86  if (hit2->isValid() && !valid) {
87  if (debug_) std::cout << " replacing old invalid hit on detid " << id2() << std::endl;
88  hit2 = hit; shared = true; break;
89  }
90  if (debug_) std::cout << " discared as additional hit on layer that already contains hit with detid " << id() << std::endl;
91  shared = true; break;
92  }
93  if (shared) continue;
94  hits.push_back(hit);
95  }
96 #endif
97 
98  math::XYZVector p = (inner.innerMomentum() + outer.outerMomentum());
99  GlobalVector v(p.x(), p.y(), p.z());
100  if (!useInnermostState_) v *= -1;
101 
103  unsigned int nhits = hits.size();
104  ownHits.reserve(nhits);
105 
106 #if TRACK_SORT
107  if (!useInnermostState_) std::reverse(hits.begin(), hits.end());
108  foreach(const TrackingRecHit * hit, hits) {
109  ownHits.push_back(*hit);
110  }
111 #elif DET_SORT
112  // OLD METHOD, sometimes fails
113  std::sort(hits.begin(), hits.end(), MomentumSort(v, &*theGeometry));
114 
115  foreach(const TrackingRecHit * hit, hits) {
116  ownHits.push_back(*hit);
117  }
118 #else
119  // NEW sort, more accurate
120  std::vector<TransientTrackingRecHit::RecHitPointer> ttrh(nhits);
121  for (unsigned int i = 0; i < nhits; ++i) ttrh[i] = theBuilder->build(hits[i]);
122  std::sort(ttrh.begin(), ttrh.end(), GlobalMomentumSort(v));
123 
124  foreach(const TransientTrackingRecHit::RecHitPointer & hit, ttrh) {
125  ownHits.push_back(*hit->hit());
126  }
127 #endif
128 
129  PTrajectoryStateOnDet state;
131  if (useInnermostState_) {
132  pdir = alongMomentum;
133  if ((inner.outerPosition()-inner.innerPosition()).Dot(inner.momentum()) >= 0) {
134  // use inner state
136  state = trajectoryStateTransform::persistentState( originalTsosIn, DetId(inner.innerDetId()) );
137  } else {
138  // use outer state
140  state = trajectoryStateTransform::persistentState( originalTsosOut, DetId(inner.outerDetId()) );
141  }
142  } else {
143  pdir = oppositeToMomentum;
144  if ((outer.outerPosition()-inner.innerPosition()).Dot(inner.momentum()) >= 0) {
145  // use outer state
147  state = trajectoryStateTransform::persistentState( originalTsosOut, DetId(outer.outerDetId()) );
148  } else {
149  // use inner state
151  state = trajectoryStateTransform::persistentState( originalTsosIn, DetId(outer.innerDetId()) );
152  }
153 
154  }
156  return TrackCandidate(ownHits, seed, state, (useInnermostState_ ? inner : outer).seedRef());
157 }
158 
159 
161 {
162  const GeomDet *det1 = geom_->idToDet(hit1->geographicalId());
163  const GeomDet *det2 = geom_->idToDet(hit2->geographicalId());
164  GlobalPoint p1 = det1->toGlobal(LocalPoint(0,0,0));
165  GlobalPoint p2 = det2->toGlobal(LocalPoint(0,0,0));
166  return (p2 - p1).dot(dir_) > 0;
167 }
168 
170 {
171  GlobalPoint p1 = hit1->isValid() ? hit1->globalPosition() : hit1->det()->position();
172  GlobalPoint p2 = hit2->isValid() ? hit2->globalPosition() : hit2->det()->position();
173  return (p2 - p1).dot(dir_) > 0;
174 }
edm::ESHandle< TrackerTopology > theTrkTopo
Definition: TrackMerger.h:27
int i
Definition: DBlmapReader.cc:9
const TrackerGeometry * geom_
Definition: TrackMerger.h:42
static const char dir_[]
std::string theBuilderName
Definition: TrackMerger.h:25
TrajectoryStateOnSurface outerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field, bool withErr=true)
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:119
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:52
size_type size() const
Definition: OwnVector.h:254
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
virtual bool sharesInput(const TrackingRecHit *other, SharedInputType what) const
PropagationDirection
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:723
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:65
void init(const edm::EventSetup &iSetup)
Definition: TrackMerger.cc:30
TrackCandidate merge(const reco::Track &inner, const reco::Track &outer) const
Definition: TrackMerger.cc:38
bool operator()(const TransientTrackingRecHit::RecHitPointer &hit1, const TransientTrackingRecHit::RecHitPointer &hit2) const
Definition: TrackMerger.cc:169
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:55
void push_back(D *&d)
Definition: OwnVector.h:280
bool operator()(const TrackingRecHit *hit1, const TrackingRecHit *hit2) const
Definition: TrackMerger.cc:160
unsigned int outerDetId() const
DetId of the detector on which surface the outermost state is located.
Definition: Track.h:94
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:104
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
double p2[4]
Definition: TauolaWrapper.h:90
bool useInnermostState_
Definition: TrackMerger.h:23
std::shared_ptr< TrackingRecHit const > RecHitPointer
Definition: DetId.h:18
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:70
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
const T & get() const
Definition: EventSetup.h:55
bool isValid() const
edm::ESHandle< TransientTrackingRecHitBuilder > theBuilder
Definition: TrackMerger.h:26
double p1[4]
Definition: TauolaWrapper.h:89
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:60
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
tuple cout
Definition: gather_cfg.py:121
DetId geographicalId() const
volatile std::atomic< bool > shutdown_flag false
bool isValid() const
Definition: ESHandle.h:47
TrackMerger(const edm::ParameterSet &iConfig)
Definition: TrackMerger.cc:19
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:99
edm::ESHandle< MagneticField > theMagField
Definition: TrackMerger.h:22
void reserve(size_t)
Definition: OwnVector.h:274
edm::ESHandle< TrackerGeometry > theGeometry
Definition: TrackMerger.h:21
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
virtual const TrackerGeomDet * idToDet(DetId) const
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field, bool withErr=true)
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:109