CMS 3D CMS Logo

storeTracks.h
Go to the documentation of this file.
1 #ifndef RecoTracker_PixelTrackFitting_plugins_storeTracks_h
2 #define RecoTracker_PixelTrackFitting_plugins_storeTracks_h
3 
7 
14 
17 
18 template <typename Ev, typename TWH>
19 void storeTracks(Ev& ev, const TWH& tracksWithHits, const TrackerTopology& ttopo) {
20  auto tracks = std::make_unique<reco::TrackCollection>();
21  auto recHits = std::make_unique<TrackingRecHitCollection>();
22  auto trackExtras = std::make_unique<reco::TrackExtraCollection>();
23 
24  int cc = 0, nTracks = tracksWithHits.size();
25 
26  trackExtras->resize(nTracks);
27  tracks->reserve(nTracks);
28  recHits->reserve(4 * nTracks);
29 
30  for (int i = 0; i < nTracks; i++) {
31  reco::Track* track = tracksWithHits[i].first;
32  const auto& hits = tracksWithHits[i].second;
33 
34  for (unsigned int k = 0; k < hits.size(); k++) {
35  auto* hit = hits[k]->clone(); // need to clone (at least if from SoA)
36  track->appendHitPattern(*hit, ttopo);
37  recHits->push_back(hit);
38  }
39  tracks->push_back(*track);
40  delete track;
41  }
42 
43  LogDebug("TrackProducer") << "put the collection of TrackingRecHit in the event"
44  << "\n";
46 
48  for (int k = 0; k < nTracks; k++) {
49  auto& aTrackExtra = (*trackExtras)[k];
50 
51  //fill the TrackExtra with TrackingRecHitRef
52  unsigned int nHits = (*tracks)[k].numberOfValidHits();
53  aTrackExtra.setHits(hitCollProd, cc, nHits);
54  cc += nHits;
55  AlgebraicVector5 v = AlgebraicVector5(0, 0, 0, 0, 0);
58  aTrackExtra.setTrajParams(std::move(trajParams), std::move(chi2s));
59  }
60 
61  LogDebug("TrackProducer") << "put the collection of TrackExtra in the event"
62  << "\n";
64 
65  for (int k = 0; k < nTracks; k++) {
66  const reco::TrackExtraRef theTrackExtraRef(ohTE, k);
67  (*tracks)[k].setExtra(theTrackExtraRef);
68  }
69 
70  ev.put(std::move(tracks));
71 }
72 
73 #endif
std::vector< unsigned char > Chi2sFive
void storeTracks(Ev &ev, const TWH &tracksWithHits, const TrackerTopology &ttopo)
Definition: storeTracks.h:19
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
std::vector< LocalTrajectoryParameters > TrajParams
ROOT::Math::SVector< double, 5 > AlgebraicVector5
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)