-
Notifications
You must be signed in to change notification settings - Fork 4.2k
/
TrackMerger.cc
202 lines (177 loc) · 8.75 KB
/
TrackMerger.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
#include "TrackMerger.h"
#include "TrackingTools/Records/interface/TransientRecHitRecord.h"
#include "Geometry/Records/interface/IdealGeometryRecord.h"
#include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
#include "Geometry/Records/interface/TrackerTopologyRcd.h"
#include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
#include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
#define TRACK_SORT 1 // just use all hits from inner track, then append from outer outside it
#define DET_SORT 0 // sort hits using global position of detector
#define HIT_SORT 0 // sort hits using global position of hit
#include "FWCore/MessageLogger/interface/MessageLogger.h"
// #define VI_DEBUG
#ifdef VI_DEBUG
#define DPRINT(x) std::cout << x << ": "
#define PRINT std::cout
#else
#define DPRINT(x) LogTrace(x)
#define PRINT LogTrace("")
#endif
TrackMerger::TrackMerger(const edm::ParameterSet &iConfig) :
useInnermostState_(iConfig.getParameter<bool>("useInnermostState")),
theBuilderName(iConfig.getParameter<std::string>("ttrhBuilderName"))
{
}
TrackMerger::~TrackMerger()
{
}
void TrackMerger::init(const edm::EventSetup &iSetup)
{
iSetup.get<TrackerDigiGeometryRecord>().get(theGeometry);
iSetup.get<IdealMagneticFieldRecord>().get(theMagField);
iSetup.get<TransientRecHitRecord>().get(theBuilderName,theBuilder);
iSetup.get<TrackerTopologyRcd>().get(theTrkTopo);
}
TrackCandidate TrackMerger::merge(const reco::Track &inner, const reco::Track &outer, DuplicateTrackType duplicateType) const
{
DPRINT("TrackMerger") << std::abs(inner.eta()) << " merging " << inner.algo() << '/' << outer.algo() << ' ' << inner.eta() << '/' << outer.eta()<< std::endl;
std::vector<const TrackingRecHit *> hits;
hits.reserve(inner.recHitsSize() + outer.recHitsSize());
DPRINT("TrackMerger") << "Inner track hits: " << std::endl;
for (auto it = inner.recHitsBegin(), ed = inner.recHitsEnd(); it != ed; ++it) {
hits.push_back(&**it);
if (debug_) {
DetId id(hits.back()->geographicalId());
PRINT << " subdet " << id.subdetId() << " layer " << theTrkTopo->layer(id) << " valid " << hits.back()->isValid() << " detid: " << id() << std::endl;
}
}
DPRINT("TrackMerger") << "Outer track hits: " << std::endl;
if(duplicateType == DuplicateTrackType::Disjoint) {
#if TRACK_SORT
DetId lastId(hits.back()->geographicalId());
int lastSubdet = lastId.subdetId();
unsigned int lastLayer = theTrkTopo->layer(lastId);
for (auto it = outer.recHitsBegin(), ed = outer.recHitsEnd(); it != ed; ++it) {
const TrackingRecHit *hit = &**it;
DetId id(hit->geographicalId());
int thisSubdet = id.subdetId();
if (thisSubdet > lastSubdet || (thisSubdet == lastSubdet && theTrkTopo->layer(id) > lastLayer)) {
hits.push_back(hit);
PRINT << " adding subdet " << id.subdetId() << " layer " << theTrkTopo->layer(id) << " valid " << hit->isValid() << " detid: " << id() << std::endl;
} else {
PRINT << " skipping subdet " << thisSubdet << " layer " << theTrkTopo->layer(id) << " valid " << hit->isValid() << " detid: " << id() << std::endl;
}
}
#else
addSecondTrackHits(hits, outer);
#endif
}
else if(duplicateType == DuplicateTrackType::Overlapping) {
addSecondTrackHits(hits, outer);
}
math::XYZVector p = (inner.innerMomentum() + outer.outerMomentum());
GlobalVector v(p.x(), p.y(), p.z());
if (!useInnermostState_) v *= -1;
TrackCandidate::RecHitContainer ownHits;
unsigned int nhits = hits.size();
ownHits.reserve(nhits);
if(duplicateType == DuplicateTrackType::Disjoint) {
#if TRACK_SORT
if (!useInnermostState_) std::reverse(hits.begin(), hits.end());
for(auto hit : hits) ownHits.push_back(*hit);
#elif DET_SORT
// OLD METHOD, sometimes fails
std::sort(hits.begin(), hits.end(), MomentumSort(v, &*theGeometry));
for(auto hit : hits) ownHits.push_back(*hit);
#else
sortByHitPosition(v, hits, ownHits);
#endif
}
else if(duplicateType == DuplicateTrackType::Overlapping) {
sortByHitPosition(v, hits, ownHits);
}
PTrajectoryStateOnDet state;
PropagationDirection pdir;
if (useInnermostState_) {
pdir = alongMomentum;
if ((inner.outerPosition()-inner.innerPosition()).Dot(inner.momentum()) >= 0) {
// use inner state
TrajectoryStateOnSurface originalTsosIn(trajectoryStateTransform::innerStateOnSurface(inner, *theGeometry, &*theMagField));
state = trajectoryStateTransform::persistentState( originalTsosIn, DetId(inner.innerDetId()) );
} else {
// use outer state
TrajectoryStateOnSurface originalTsosOut(trajectoryStateTransform::outerStateOnSurface(inner, *theGeometry, &*theMagField));
state = trajectoryStateTransform::persistentState( originalTsosOut, DetId(inner.outerDetId()) );
}
} else {
pdir = oppositeToMomentum;
if ((outer.outerPosition()-inner.innerPosition()).Dot(inner.momentum()) >= 0) {
// use outer state
TrajectoryStateOnSurface originalTsosOut(trajectoryStateTransform::outerStateOnSurface(outer, *theGeometry, &*theMagField));
state = trajectoryStateTransform::persistentState( originalTsosOut, DetId(outer.outerDetId()) );
} else {
// use inner state
TrajectoryStateOnSurface originalTsosIn(trajectoryStateTransform::innerStateOnSurface(outer, *theGeometry, &*theMagField));
state = trajectoryStateTransform::persistentState( originalTsosIn, DetId(outer.innerDetId()) );
}
}
TrajectorySeed seed(state, TrackCandidate::RecHitContainer(), pdir);
TrackCandidate ret(ownHits, seed, state, (useInnermostState_ ? inner : outer).seedRef());
ret.setStopReason((uint8_t)(useInnermostState_ ? inner : outer).stopReason());
return ret;
}
void TrackMerger::addSecondTrackHits(std::vector<const TrackingRecHit *>& hits,
const reco::Track& outer) const {
size_t nHitsFirstTrack = hits.size();
for (auto it = outer.recHitsBegin(), ed = outer.recHitsEnd(); it != ed; ++it) {
const TrackingRecHit *hit = &**it;
DetId id(hit->geographicalId());
const auto lay = theTrkTopo->layer(id);
bool shared = false;
bool valid = hit->isValid();
PRINT << " subdet " << id.subdetId() << " layer " << theTrkTopo->layer(id) << " valid " << valid << " detid: " << id() << std::endl;
size_t iHit = 0;
for ( auto hit2 : hits) {
++iHit; if (iHit > nHitsFirstTrack) break;
DetId id2 = hit2->geographicalId();
if (id.subdetId() != id2.subdetId()) continue;
if (theTrkTopo->layer(id2) != lay) continue;
if (hit->sharesInput(hit2, TrackingRecHit::all)) {
PRINT << " discared as duplicate of other hit" << id() << std::endl;
shared = true; break;
}
if (hit2->isValid() && !valid) {
PRINT << " replacing old invalid hit on detid " << id2() << std::endl;
hit2 = hit; shared = true; break;
}
PRINT << " discared as additional hit on layer that already contains hit with detid " << id() << std::endl;
shared = true; break;
}
if (shared) continue;
hits.push_back(hit);
}
}
void TrackMerger::sortByHitPosition(const GlobalVector& v,
const std::vector<const TrackingRecHit *>& hits,
TrackCandidate::RecHitContainer& ownHits) const {
// NEW sort, more accurate
unsigned int nhits = hits.size();
std::vector<TransientTrackingRecHit::RecHitPointer> ttrh(nhits);
for (unsigned int i = 0; i < nhits; ++i) ttrh[i] = theBuilder->build(hits[i]);
std::sort(ttrh.begin(), ttrh.end(), GlobalMomentumSort(v));
for(auto hit : ttrh) ownHits.push_back(*hit);
}
bool TrackMerger::MomentumSort::operator()(const TrackingRecHit *hit1, const TrackingRecHit *hit2) const
{
const GeomDet *det1 = geom_->idToDet(hit1->geographicalId());
const GeomDet *det2 = geom_->idToDet(hit2->geographicalId());
GlobalPoint p1 = det1->toGlobal(LocalPoint(0,0,0));
GlobalPoint p2 = det2->toGlobal(LocalPoint(0,0,0));
return (p2 - p1).dot(dir_) > 0;
}
bool TrackMerger::GlobalMomentumSort::operator()(const TransientTrackingRecHit::RecHitPointer &hit1, const TransientTrackingRecHit::RecHitPointer &hit2) const
{
GlobalPoint p1 = hit1->isValid() ? hit1->globalPosition() : hit1->det()->position();
GlobalPoint p2 = hit2->isValid() ? hit2->globalPosition() : hit2->det()->position();
return (p2 - p1).dot(dir_) > 0;
}