CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelTrackProducer.cc
Go to the documentation of this file.
1 #include "PixelTrackProducer.h"
2 
8 
13 
16 
17 #include <vector>
18 
19 using namespace pixeltrackfitting;
20 using edm::ParameterSet;
21 
23  : theReconstruction(cfg, consumesCollector())
24 {
25  edm::LogInfo("PixelTrackProducer")<<" construction...";
26  produces<reco::TrackCollection>();
27  produces<TrackingRecHitCollection>();
28  produces<reco::TrackExtraCollection>();
29 }
30 
32 
34 {
36 }
37 
39 {
41 }
42 
44 {
45  LogDebug("PixelTrackProducer, produce")<<"event# :"<<ev.id();
46 
48  theReconstruction.run(tracks,ev,es);
49 
51  es.get<TrackerTopologyRcd>().get(httopo);
52 
53  // store tracks
54  store(ev, tracks, *httopo);
55 }
56 
57 void PixelTrackProducer::store(edm::Event& ev, const TracksWithTTRHs& tracksWithHits, const TrackerTopology& ttopo)
58 {
59  std::auto_ptr<reco::TrackCollection> tracks(new reco::TrackCollection());
60  std::auto_ptr<TrackingRecHitCollection> recHits(new TrackingRecHitCollection());
61  std::auto_ptr<reco::TrackExtraCollection> trackExtras(new reco::TrackExtraCollection());
62 
63  int cc = 0, nTracks = tracksWithHits.size();
64 
65  for (int i = 0; i < nTracks; i++)
66  {
67  reco::Track* track = tracksWithHits.at(i).first;
68  const SeedingHitSet& hits = tracksWithHits.at(i).second;
69 
70  for (unsigned int k = 0; k < hits.size(); k++)
71  {
72  TrackingRecHit *hit = hits[k]->hit()->clone();
73 
74  track->appendHitPattern(*hit, ttopo);
75  recHits->push_back(hit);
76  }
77  tracks->push_back(*track);
78  delete track;
79 
80  }
81 
82  LogDebug("TrackProducer") << "put the collection of TrackingRecHit in the event" << "\n";
84 
86  for (int k = 0; k < nTracks; k++)
87  {
88  reco::TrackExtra theTrackExtra{};
89 
90  //fill the TrackExtra with TrackingRecHitRef
91  unsigned int nHits = tracks->at(k).numberOfValidHits();
92  theTrackExtra.setHits(hitCollProd, cc, nHits);
93  cc +=nHits;
94  trackExtras->push_back(theTrackExtra);
95  }
96 
97  LogDebug("TrackProducer") << "put the collection of TrackExtra in the event" << "\n";
99 
100  for (int k = 0; k < nTracks; k++)
101  {
102  const reco::TrackExtraRef theTrackExtraRef(ohTE,k);
103  (tracks->at(k)).setExtra(theTrackExtraRef);
104  }
105 
106  ev.put(tracks);
107 
108 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
const unsigned int nTracks(const reco::Vertex &sv)
void init(const edm::EventSetup &es)
tuple cfg
Definition: looper.py:293
void setHits(TrackingRecHitRefProd const &prod, unsigned firstH, unsigned int nH)
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
bool ev
virtual void endRun(const edm::Run &run, const edm::EventSetup &es) override
std::vector< TrackWithTTRHs > TracksWithTTRHs
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
virtual void beginRun(const edm::Run &run, const edm::EventSetup &es) override
virtual void produce(edm::Event &ev, const edm::EventSetup &es) override
void run(pixeltrackfitting::TracksWithTTRHs &tah, edm::Event &ev, const edm::EventSetup &es)
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:11
edm::OwnVector< TrackingRecHit > TrackingRecHitCollection
collection of TrackingRecHits
tuple tracks
Definition: testEve_cfg.py:39
const T & get() const
Definition: EventSetup.h:56
bool appendHitPattern(const TrackingRecHit &hit, const TrackerTopology &ttopo)
append a single hit to the HitPattern
Definition: TrackBase.h:455
PixelTrackReconstruction theReconstruction
edm::EventID id() const
Definition: EventBase.h:59
unsigned int size() const
Definition: SeedingHitSet.h:44
void store(edm::Event &ev, const pixeltrackfitting::TracksWithTTRHs &selectedTracks, const TrackerTopology &ttopo)
PixelTrackProducer(const edm::ParameterSet &conf)
Definition: Run.h:43