CMS 3D CMS Logo

TransientTrackBuilder.cc
Go to the documentation of this file.
6 
7 
8 using namespace reco;
9 using namespace std;
10 using namespace edm;
11 
12 namespace {
13  constexpr float fakeBeamSpotTimeWidth = 0.175f;
14 }
15 
17  return TransientTrack(*t, theField, theTrackingGeometry);
18 }
19 
21  return TransientTrack(t, theField, theTrackingGeometry);
22 }
23 
25  return TransientTrack(new GsfTransientTrack(*t, theField, theTrackingGeometry));
26 }
27 
29  return TransientTrack(new GsfTransientTrack(t, theField, theTrackingGeometry));
30 }
31 
33  return TransientTrack(*t, theField, theTrackingGeometry);
34 }
35 
37  return TransientTrack(t, theField, theTrackingGeometry);
38 }
39 
41  return TransientTrack(*t, theField, theTrackingGeometry);
42 }
43 
45  return TransientTrack(t, theField, theTrackingGeometry);
46 }
47 
48 
50  return TransientTrack(new GsfTransientTrack(*t, theField, theTrackingGeometry));
51 }
52 
54  return TransientTrack(new GsfTransientTrack(t, theField, theTrackingGeometry));
55 }
56 
57 vector<TransientTrack>
59 {
60  vector<TransientTrack> ttVect;
61  ttVect.reserve((*trkColl).size());
62  for (unsigned int i = 0; i < (*trkColl).size() ; i++) {
63  ttVect.push_back(TransientTrack(TrackRef(trkColl, i), theField, theTrackingGeometry));
64  }
65  return ttVect;
66 }
67 
68 vector<TransientTrack>
70 {
71  vector<TransientTrack> ttVect;
72  ttVect.reserve((*trkColl).size());
73  for (unsigned int i = 0; i < (*trkColl).size() ; i++) {
74  ttVect.push_back( TransientTrack(
75  new GsfTransientTrack(GsfTrackRef(trkColl, i), theField, theTrackingGeometry)) );
76  }
77  return ttVect;
78 }
79 
80 vector<TransientTrack>
82 {
83  vector<TransientTrack> ttVect;
84  ttVect.reserve((*trkColl).size());
85  for (unsigned int i = 0; i < (*trkColl).size() ; i++) {
86  const Track * trk = &(*trkColl)[i];
87  const GsfTrack * gsfTrack = dynamic_cast<const GsfTrack *>(trk);
88  if (gsfTrack) {
89  ttVect.push_back( TransientTrack(
90  new GsfTransientTrack(RefToBase<Track>(trkColl, i).castTo<GsfTrackRef>(), theField, theTrackingGeometry)) );
91  } else { // gsf
92  ttVect.push_back(TransientTrack(RefToBase<Track>(trkColl, i).castTo<TrackRef>(), theField, theTrackingGeometry));
93  }
94  }
95  return ttVect;
96 }
97 
98 vector<TransientTrack>
100  const edm::ValueMap<float>& trackTimes,
101  const edm::ValueMap<float>& trackTimeResos ) const
102 {
103  vector<TransientTrack> ttVect;
104  ttVect.reserve((*trkColl).size());
105  for (unsigned int i = 0; i < (*trkColl).size() ; i++) {
106  TrackRef ref(trkColl, i);
107  double time = trackTimes[ref];
108  double timeReso = trackTimeResos[ref];
109  timeReso = ( timeReso > 1e-6 ? timeReso : fakeBeamSpotTimeWidth ); // make the error much larger than the BS time width
110  if( edm::isNotFinite(time) ) {
111  time = 0.0;
112  timeReso = 2.0*fakeBeamSpotTimeWidth;
113  }
114  ttVect.push_back(TransientTrack(ref, time, timeReso, theField, theTrackingGeometry));
115  }
116  return ttVect;
117 }
118 
119 vector<TransientTrack>
121  const edm::ValueMap<float>& trackTimes,
122  const edm::ValueMap<float>& trackTimeResos ) const
123 {
124  vector<TransientTrack> ttVect;
125  ttVect.reserve((*trkColl).size());
126  for (unsigned int i = 0; i < (*trkColl).size() ; i++) {
127  GsfTrackRef ref(trkColl, i);
128  double time = trackTimes[ref];
129  double timeReso = trackTimeResos[ref];
130  timeReso = ( timeReso > 1e-6 ? timeReso : fakeBeamSpotTimeWidth ); // make the error much larger than the BS time width
131  if( edm::isNotFinite(time) ) {
132  time = 0.0;
133  timeReso = 2.0*fakeBeamSpotTimeWidth;
134  }
135  ttVect.push_back( TransientTrack(
136  new GsfTransientTrack(ref, time, timeReso, theField, theTrackingGeometry)) );
137  }
138  return ttVect;
139 }
140 
141 vector<TransientTrack>
143  const edm::ValueMap<float>& trackTimes,
144  const edm::ValueMap<float>& trackTimeResos ) const
145 {
146  vector<TransientTrack> ttVect;
147  ttVect.reserve((*trkColl).size());
148  for (unsigned int i = 0; i < (*trkColl).size() ; i++) {
149  const Track * trk = &(*trkColl)[i];
150  const GsfTrack * gsfTrack = dynamic_cast<const GsfTrack *>(trk);
151  if (gsfTrack) {
152  GsfTrackRef ref = RefToBase<Track>(trkColl, i).castTo<GsfTrackRef>();
153  double time = trackTimes[ref];
154  double timeReso = trackTimeResos[ref];
155  timeReso = ( timeReso > 1e-6 ? timeReso : fakeBeamSpotTimeWidth ); // make the error much larger than the BS time width
156  if( edm::isNotFinite(time) ) {
157  time = 0.0;
158  timeReso = 2.0*fakeBeamSpotTimeWidth;
159  }
160  ttVect.push_back( TransientTrack(
161  new GsfTransientTrack(RefToBase<Track>(trkColl, i).castTo<GsfTrackRef>(), time, timeReso, theField, theTrackingGeometry)) );
162  } else { // gsf
163  TrackRef ref = RefToBase<Track>(trkColl, i).castTo<TrackRef>();
164  double time = trackTimes[ref];
165  double timeReso = trackTimeResos[ref];
166  timeReso = ( timeReso > 1e-6 ? timeReso : fakeBeamSpotTimeWidth ); // make the error much larger than the BS time width
167  if( edm::isNotFinite(time) ) {
168  time = 0.0;
169  timeReso = 2.0*fakeBeamSpotTimeWidth;
170  }
171  ttVect.push_back(TransientTrack(RefToBase<Track>(trkColl, i).castTo<TrackRef>(), time, timeReso, theField, theTrackingGeometry));
172  }
173  }
174  return ttVect;
175 }
176 
177 vector<TransientTrack>
179  const reco::BeamSpot & beamSpot) const
180 {
181  vector<TransientTrack> ttVect = build(trkColl);
182  for (unsigned int i = 0; i < ttVect.size() ; i++) {
183  ttVect[i].setBeamSpot(beamSpot);
184  }
185  return ttVect;
186 }
187 
188 vector<TransientTrack>
190  const reco::BeamSpot & beamSpot) const
191 {
192  vector<TransientTrack> ttVect = build(trkColl);
193  for (unsigned int i = 0; i < ttVect.size() ; i++) {
194  ttVect[i].setBeamSpot(beamSpot);
195  }
196  return ttVect;
197 }
198 
199 vector<TransientTrack>
201  const reco::BeamSpot & beamSpot) const
202 {
203  vector<TransientTrack> ttVect = build(trkColl);
204  for (unsigned int i = 0; i < ttVect.size() ; i++) {
205  ttVect[i].setBeamSpot(beamSpot);
206  }
207  return ttVect;
208 }
209 
210 vector<TransientTrack>
212  const reco::BeamSpot & beamSpot,
213  const edm::ValueMap<float>& trackTimes,
214  const edm::ValueMap<float>& trackTimeResos ) const
215 {
216  vector<TransientTrack> ttVect = build(trkColl, trackTimes, trackTimeResos );
217  for (unsigned int i = 0; i < ttVect.size() ; i++) {
218  ttVect[i].setBeamSpot(beamSpot);
219  }
220  return ttVect;
221 }
222 
223 vector<TransientTrack>
225  const reco::BeamSpot & beamSpot,
226  const edm::ValueMap<float>& trackTimes,
227  const edm::ValueMap<float>& trackTimeResos ) const
228 {
229  vector<TransientTrack> ttVect = build(trkColl, trackTimes, trackTimeResos);
230  for (unsigned int i = 0; i < ttVect.size() ; i++) {
231  ttVect[i].setBeamSpot(beamSpot);
232  }
233  return ttVect;
234 }
235 
236 vector<TransientTrack>
238  const reco::BeamSpot & beamSpot,
239  const edm::ValueMap<float>& trackTimes,
240  const edm::ValueMap<float>& trackTimeResos ) const
241 {
242  vector<TransientTrack> ttVect = build(trkColl, trackTimes, trackTimeResos);
243  for (unsigned int i = 0; i < ttVect.size() ; i++) {
244  ttVect[i].setBeamSpot(beamSpot);
245  }
246  return ttVect;
247 }
248 
250  return TransientTrack(new TransientTrackFromFTS(fts));
251 }
edm::Ref< GsfTrackCollection > GsfTrackRef
persistent reference to a GsfTrack
Definition: GsfTrackFwd.h:13
reco::TransientTrack build(const reco::Track *p) const
#define constexpr
bool isNotFinite(T x)
Definition: isFinite.h:10
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
fixed size matrix
HLT enums.