CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelTracksProducer.cc
Go to the documentation of this file.
1 #include <memory>
3 
7 
10 
12 
13 //Pixel Specific stuff
16 
19 
24 
29 #include <vector>
30 
31 using namespace pixeltrackfitting;
32 
34  theFitter(0),
35  theFilter(0),
36  theRegionProducer(0)
37 {
38 
39  produces<reco::TrackCollection>();
40  produces<TrackingRecHitCollection>();
41  produces<reco::TrackExtraCollection>();
42 
43  const edm::ParameterSet& regfactoryPSet = conf.getParameter<edm::ParameterSet>("RegionFactoryPSet");
44  std::string regfactoryName = regfactoryPSet.getParameter<std::string>("ComponentName");
45  theRegionProducer = TrackingRegionProducerFactory::get()->create(regfactoryName,
46  regfactoryPSet, consumesCollector());
47 
48  const edm::ParameterSet& fitterPSet = conf.getParameter<edm::ParameterSet>("FitterPSet");
49  std::string fitterName = fitterPSet.getParameter<std::string>("ComponentName");
50  theFitter = PixelFitterFactory::get()->create( fitterName, fitterPSet);
51 
53  const edm::ParameterSet& filterPSet = conf.getParameter<edm::ParameterSet>("FilterPSet");
54  std::string filterName = filterPSet.getParameter<std::string>("ComponentName");
55  theFilter = PixelTrackFilterFactory::get()->create( filterName, filterPSet, iC);
56 
57  // The name of the seed producer
58  seedProducer = conf.getParameter<edm::InputTag>("SeedProducer");
59  seedProducerToken = consumes<TrajectorySeedCollection>(seedProducer);
60 
61 }
62 
63 
64 // Virtual destructor needed.
66 
67  delete theFilter;
68  delete theFitter;
69  delete theRegionProducer;
70 
71 }
72 
73 
74 // Functions that gets called by framework every event
75 void
77 
78  std::auto_ptr<reco::TrackCollection> tracks(new reco::TrackCollection);
79  std::auto_ptr<TrackingRecHitCollection> recHits(new TrackingRecHitCollection);
80  std::auto_ptr<reco::TrackExtraCollection> trackExtras(new reco::TrackExtraCollection);
81  typedef std::vector<const TrackingRecHit *> RecHits;
82 
83  TracksWithRecHits pixeltracks;
84  TracksWithRecHits cleanedTracks;
85 
87  e.getByToken(seedProducerToken,theSeeds);
88 
89  // No seed -> output an empty track collection
90  if(theSeeds->size() == 0) {
91  e.put(tracks);
92  e.put(recHits);
93  e.put(trackExtras);
94  return;
95  }
96 
97  //only one region Global, but it is called at every event...
98  //maybe there is a smarter way to set it only once
99  //NEED TO FIX
100  typedef std::vector<TrackingRegion* > Regions;
101  typedef Regions::const_iterator IR;
102  Regions regions = theRegionProducer->regions(e,es);
103  for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) {
104  const TrackingRegion & region = **ir;
105 
106  // Loop over the seeds
107  TrajectorySeedCollection::const_iterator aSeed = theSeeds->begin();
108  TrajectorySeedCollection::const_iterator lastSeed = theSeeds->end();
109  for ( ; aSeed!=lastSeed; ++aSeed ) {
110 
111  // Find the first hit and last hit of the Seed
112  TrajectorySeed::range theSeedingRecHitRange = aSeed->recHits();
113  edm::OwnVector<TrackingRecHit>::const_iterator aSeedingRecHit = theSeedingRecHitRange.first;
114  edm::OwnVector<TrackingRecHit>::const_iterator theLastSeedingRecHit = theSeedingRecHitRange.second;
115 
116  // Loop over the rechits
117  std::vector<const TrackingRecHit*> TripletHits(3,static_cast<const TrackingRecHit*>(0));
118  for ( unsigned i=0; aSeedingRecHit!=theLastSeedingRecHit; ++i,++aSeedingRecHit )
119  TripletHits[i] = &(*aSeedingRecHit);
120 
121  // fitting the triplet
122  reco::Track* track = theFitter->run(es, TripletHits, region);
123 
124  // decide if track should be skipped according to filter
125  if ( ! (*theFilter)(track) ) {
126  delete track;
127  continue;
128  }
129 
130  // add tracks
131  pixeltracks.push_back(TrackWithRecHits(track, TripletHits));
132 
133  }
134  }
135 
136  int cc=0;
137  int nTracks = pixeltracks.size();
138  for (int i = 0; i < nTracks; ++i) {
139 
140  reco::Track* track = pixeltracks.at(i).first;
141  const RecHits & hits = pixeltracks.at(i).second;
142 
143  for (unsigned int k = 0; k < hits.size(); k++) {
144  TrackingRecHit *hit = (hits.at(k))->clone();
145  track->appendHitPattern(*hit);
146  recHits->push_back(hit);
147  }
148 
149  tracks->push_back(*track);
150  delete track;
151 
152  }
153 
156 
157  for (int k = 0; k < nTracks; ++k) {
158 
159  // reco::TrackExtra* theTrackExtra = new reco::TrackExtra();
160  reco::TrackExtra theTrackExtra;
161 
162  //fill the TrackExtra with TrackingRecHitRef
163  // unsigned int nHits = tracks->at(k).numberOfValidHits();
164  const unsigned nHits = 3; // We are dealing with triplets!
165  theTrackExtra.setHits( ohRHProd, cc, nHits);
166  cc += nHits;
167 
168  trackExtras->push_back(theTrackExtra);
169  //trackExtras->push_back(*theTrackExtra);
170  //delete theTrackExtra;
171  }
172 
174 
175  for (int k = 0; k < nTracks; k++) {
176 
177  const reco::TrackExtraRef theTrackExtraRef(ohTE,k);
178  (tracks->at(k)).setExtra(theTrackExtraRef);
179 
180  }
181 
182  e.put(tracks);
183 
184  // Avoid a memory leak !
185  unsigned nRegions = regions.size();
186  for ( unsigned iRegions=0; iRegions<nRegions; ++iRegions ) {
187  delete regions[iRegions];
188  }
189 
190 }
191 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
edm::InputTag seedProducer
void setHits(TrackingRecHitRefProd const &prod, unsigned firstH, unsigned int nH)
PixelTracksProducer(const edm::ParameterSet &conf)
virtual reco::Track * run(const edm::EventSetup &es, const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region) const
Definition: PixelFitter.h:17
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
virtual void produce(edm::Event &ev, const edm::EventSetup &es)
TrackingRegionProducer * theRegionProducer
edm::EDGetTokenT< TrajectorySeedCollection > seedProducerToken
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
const PixelTrackFilter * theFilter
std::pair< const_iterator, const_iterator > range
bool appendHitPattern(const TrackingRecHit &hit)
append a single hit to the HitPattern
Definition: TrackBase.h:421
std::pair< reco::Track *, std::vector< const TrackingRecHit * > > TrackWithRecHits
tuple conf
Definition: dbtoconf.py:185
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:11
tuple tracks
Definition: testEve_cfg.py:39
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135
const PixelFitter * theFitter
std::vector< TrackWithRecHits > TracksWithRecHits
T get(const Candidate &c)
Definition: component.h:55
virtual std::vector< TrackingRegion * > regions(const edm::Event &ev, const edm::EventSetup &es) const =0