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 
7 
12 
13 #include <vector>
14 
15 using namespace pixeltrackfitting;
16 using edm::ParameterSet;
17 
19  : theReconstruction(cfg, consumesCollector())
20 {
21  edm::LogInfo("PixelTrackProducer")<<" construction...";
22  produces<reco::TrackCollection>();
23  produces<TrackingRecHitCollection>();
24  produces<reco::TrackExtraCollection>();
25 }
26 
28 
30 {
32 }
33 
35 {
37 }
38 
40 {
41  LogDebug("PixelTrackProducer, produce")<<"event# :"<<ev.id();
42 
44  theReconstruction.run(tracks,ev,es);
45 
46  // store tracks
47  store(ev, tracks);
48 }
49 
50 void PixelTrackProducer::store(edm::Event& ev, const TracksWithTTRHs& tracksWithHits)
51 {
52  std::auto_ptr<reco::TrackCollection> tracks(new reco::TrackCollection());
53  std::auto_ptr<TrackingRecHitCollection> recHits(new TrackingRecHitCollection());
54  std::auto_ptr<reco::TrackExtraCollection> trackExtras(new reco::TrackExtraCollection());
55 
56  int cc = 0, nTracks = tracksWithHits.size();
57 
58  for (int i = 0; i < nTracks; i++)
59  {
60  reco::Track* track = tracksWithHits.at(i).first;
61  const SeedingHitSet& hits = tracksWithHits.at(i).second;
62 
63  for (unsigned int k = 0; k < hits.size(); k++)
64  {
65  TrackingRecHit *hit = hits[k]->hit()->clone();
66 
67  track->setHitPattern(*hit, k);
68  recHits->push_back(hit);
69  }
70  tracks->push_back(*track);
71  delete track;
72 
73  }
74 
75  LogDebug("TrackProducer") << "put the collection of TrackingRecHit in the event" << "\n";
77 
78 
79  for (int k = 0; k < nTracks; k++)
80  {
81  reco::TrackExtra* theTrackExtra = new reco::TrackExtra();
82 
83  //fill the TrackExtra with TrackingRecHitRef
84  unsigned int nHits = tracks->at(k).numberOfValidHits();
85  for(unsigned int i = 0; i < nHits; ++i) {
86  theTrackExtra->add(TrackingRecHitRef(ohRH,cc));
87  cc++;
88  }
89 
90  trackExtras->push_back(*theTrackExtra);
91  delete theTrackExtra;
92  }
93 
94  LogDebug("TrackProducer") << "put the collection of TrackExtra in the event" << "\n";
96 
97  for (int k = 0; k < nTracks; k++)
98  {
99  const reco::TrackExtraRef theTrackExtraRef(ohTE,k);
100  (tracks->at(k)).setExtra(theTrackExtraRef);
101  }
102 
103  ev.put(tracks);
104 
105 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
void init(const edm::EventSetup &es)
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
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:116
edm::Ref< TrackingRecHitCollection > TrackingRecHitRef
persistent reference to a TrackingRecHit
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)
int k[5][pyjets_maxn]
void setHitPattern(const C &c)
set hit patterns from vector of hit references
Definition: TrackBase.h:244
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:9
edm::OwnVector< TrackingRecHit > TrackingRecHitCollection
collection of TrackingRecHits
tuple tracks
Definition: testEve_cfg.py:39
void add(const TrackingRecHitRef &r)
add a reference to a RecHit
PixelTrackReconstruction theReconstruction
edm::EventID id() const
Definition: EventBase.h:56
unsigned int size() const
Definition: SeedingHitSet.h:41
void store(edm::Event &ev, const pixeltrackfitting::TracksWithTTRHs &selectedTracks)
PixelTrackProducer(const edm::ParameterSet &conf)
Definition: Run.h:41