CMS 3D CMS Logo

PixelTrackProducer.cc
Go to the documentation of this file.
1 #include "PixelTrackProducer.h"
2 
10 
16 
19 
20 #include <vector>
21 
22 using namespace pixeltrackfitting;
23 using edm::ParameterSet;
24 
25 PixelTrackProducer::PixelTrackProducer(const ParameterSet& cfg) : theReconstruction(cfg, consumesCollector()) {
26  edm::LogInfo("PixelTrackProducer") << " construction...";
27  produces<reco::TrackCollection>();
28  produces<TrackingRecHitCollection>();
29  produces<reco::TrackExtraCollection>();
30 }
31 
33 
36 
37  desc.add<std::string>("passLabel", "pixelTracks"); // What is this? It is not used anywhere in this code.
39 
40  descriptions.add("pixelTracks", desc);
41 }
42 
44  LogDebug("PixelTrackProducer, produce") << "event# :" << ev.id();
45 
47  theReconstruction.run(tracks, ev, es);
48 
50  es.get<TrackerTopologyRcd>().get(httopo);
51 
52  // store tracks
53  store(ev, tracks, *httopo);
54 }
55 
56 void PixelTrackProducer::store(edm::Event& ev, const TracksWithTTRHs& tracksWithHits, const TrackerTopology& ttopo) {
57  auto tracks = std::make_unique<reco::TrackCollection>();
58  auto recHits = std::make_unique<TrackingRecHitCollection>();
59  auto trackExtras = std::make_unique<reco::TrackExtraCollection>();
60 
61  int cc = 0, nTracks = tracksWithHits.size();
62 
63  for (int i = 0; i < nTracks; i++) {
64  reco::Track* track = tracksWithHits.at(i).first;
65  const SeedingHitSet& hits = tracksWithHits.at(i).second;
66 
67  for (unsigned int k = 0; k < hits.size(); k++) {
68  TrackingRecHit* hit = hits[k]->hit()->clone();
69 
70  track->appendHitPattern(*hit, ttopo);
71  recHits->push_back(hit);
72  }
73  tracks->push_back(*track);
74  delete track;
75  }
76 
77  LogDebug("TrackProducer") << "put the collection of TrackingRecHit in the event"
78  << "\n";
80 
82  for (int k = 0; k < nTracks; k++) {
83  reco::TrackExtra theTrackExtra{};
84 
85  //fill the TrackExtra with TrackingRecHitRef
86  unsigned int nHits = tracks->at(k).numberOfValidHits();
87  theTrackExtra.setHits(hitCollProd, cc, nHits);
88  cc += nHits;
89  AlgebraicVector5 v = AlgebraicVector5(0, 0, 0, 0, 0);
91  reco::TrackExtra::Chi2sFive chi2s(nHits, 0);
92  theTrackExtra.setTrajParams(std::move(trajParams), std::move(chi2s));
93  trackExtras->push_back(theTrackExtra);
94  }
95 
96  LogDebug("TrackProducer") << "put the collection of TrackExtra in the event"
97  << "\n";
99 
100  for (int k = 0; k < nTracks; k++) {
101  const reco::TrackExtraRef theTrackExtraRef(ohTE, k);
102  (tracks->at(k)).setExtra(theTrackExtraRef);
103  }
104 
105  ev.put(std::move(tracks));
106 }
#define LogDebug(id)
const unsigned int nTracks(const reco::Vertex &sv)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
std::vector< unsigned char > Chi2sFive
bool ev
std::vector< TrackWithTTRHs > TracksWithTTRHs
std::vector< LocalTrajectoryParameters > TrajParams
ROOT::Math::SVector< double, 5 > AlgebraicVector5
void produce(edm::Event &ev, const edm::EventSetup &es) override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void run(pixeltrackfitting::TracksWithTTRHs &tah, edm::Event &ev, const edm::EventSetup &es)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool appendHitPattern(const TrackingRecHit &hit, const TrackerTopology &ttopo)
append a single hit to the HitPattern
Definition: TrackBase.h:489
PixelTrackReconstruction theReconstruction
edm::EventID id() const
Definition: EventBase.h:59
T get() const
Definition: EventSetup.h:73
unsigned int size() const
Definition: SeedingHitSet.h:41
void store(edm::Event &ev, const pixeltrackfitting::TracksWithTTRHs &selectedTracks, const TrackerTopology &ttopo)
static void fillDescriptions(edm::ParameterSetDescription &desc)
~PixelTrackProducer() override
def move(src, dest)
Definition: eostools.py:511
PixelTrackProducer(const edm::ParameterSet &conf)