CMS 3D CMS Logo

PixelTracksProducer.cc
Go to the documentation of this file.
1 #include <memory>
3 
8 
11 
13 
16 
17 //Pixel Specific stuff
20 
22 
26 
31 #include <vector>
32 
33 using namespace pixeltrackfitting;
34 
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  fitterToken = consumes<PixelFitter>(conf.getParameter<edm::InputTag>("Fitter"));
49  filterToken = consumes<PixelTrackFilter>(conf.getParameter<edm::InputTag>("Filter"));
50 
51  // The name of the seed producer
52  auto seedProducer = conf.getParameter<edm::InputTag>("SeedProducer");
53  seedProducerToken = consumes<TrajectorySeedCollection>(seedProducer);
54 
55 }
56 
57 
58 // Virtual destructor needed.
60 
61  delete theRegionProducer;
62 
63 }
64 
65 
66 // Functions that gets called by framework every event
67 void
69 
70  std::unique_ptr<reco::TrackCollection> tracks(new reco::TrackCollection);
71  std::unique_ptr<TrackingRecHitCollection> recHits(new TrackingRecHitCollection);
72  std::unique_ptr<reco::TrackExtraCollection> trackExtras(new reco::TrackExtraCollection);
73  typedef std::vector<const TrackingRecHit *> RecHits;
74 
75  TracksWithRecHits pixeltracks;
76  TracksWithRecHits cleanedTracks;
77 
79  es.get<TrackerTopologyRcd>().get(httopo);
80  const TrackerTopology& ttopo = *httopo;
81 
83  e.getByToken(fitterToken, hfitter);
84  const PixelFitter& fitter = *hfitter;
85 
87  e.getByToken(filterToken, hfilter);
88  const PixelTrackFilter& theFilter = *hfilter;
89 
91  e.getByToken(seedProducerToken,theSeeds);
92 
93  // No seed -> output an empty track collection
94  if(theSeeds->size() == 0) {
95  e.put(std::move(tracks));
96  e.put(std::move(recHits));
97  e.put(std::move(trackExtras));
98  return;
99  }
100 
101  //only one region Global, but it is called at every event...
102  //maybe there is a smarter way to set it only once
103  //NEED TO FIX
104  typedef std::vector<std::unique_ptr<TrackingRegion> > Regions;
105  typedef Regions::const_iterator IR;
106  Regions regions = theRegionProducer->regions(e,es);
107  for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) {
108  const TrackingRegion & region = **ir;
109 
110  // Loop over the seeds
111  TrajectorySeedCollection::const_iterator aSeed = theSeeds->begin();
112  TrajectorySeedCollection::const_iterator lastSeed = theSeeds->end();
113  for ( ; aSeed!=lastSeed; ++aSeed ) {
114 
115  // Find the first hit and last hit of the Seed
116  TrajectorySeed::range theSeedingRecHitRange = aSeed->recHits();
117  edm::OwnVector<TrackingRecHit>::const_iterator aSeedingRecHit = theSeedingRecHitRange.first;
118  edm::OwnVector<TrackingRecHit>::const_iterator theLastSeedingRecHit = theSeedingRecHitRange.second;
119 
120  // Loop over the rechits
121  std::vector<const TrackingRecHit*> TripletHits(3,static_cast<const TrackingRecHit*>(0));
122  for ( unsigned i=0; aSeedingRecHit!=theLastSeedingRecHit; ++i,++aSeedingRecHit )
123  TripletHits[i] = &(*aSeedingRecHit);
124 
125  // fitting the triplet
126  std::unique_ptr<reco::Track> track = fitter.run(TripletHits, region);
127 
128  // decide if track should be skipped according to filter
129  if ( ! theFilter(track.get(), TripletHits) ) {
130  continue;
131  }
132 
133  // add tracks
134  pixeltracks.push_back(TrackWithRecHits(track.release(), TripletHits));
135 
136  }
137  }
138 
139  int cc=0;
140  int nTracks = pixeltracks.size();
141  for (int i = 0; i < nTracks; ++i) {
142 
143  reco::Track* track = pixeltracks.at(i).first;
144  const RecHits & hits = pixeltracks.at(i).second;
145 
146  for (unsigned int k = 0; k < hits.size(); k++) {
147  TrackingRecHit *hit = (hits.at(k))->clone();
148  track->appendHitPattern(*hit, ttopo);
149  recHits->push_back(hit);
150  }
151 
152  tracks->push_back(*track);
153  delete track;
154 
155  }
156 
159 
160  for (int k = 0; k < nTracks; ++k) {
161 
162  // reco::TrackExtra* theTrackExtra = new reco::TrackExtra();
163  reco::TrackExtra theTrackExtra;
164 
165  //fill the TrackExtra with TrackingRecHitRef
166  // unsigned int nHits = tracks->at(k).numberOfValidHits();
167  const unsigned nHits = 3; // We are dealing with triplets!
168  theTrackExtra.setHits( ohRHProd, cc, nHits);
169  cc += nHits;
170 
171  trackExtras->push_back(theTrackExtra);
172  //trackExtras->push_back(*theTrackExtra);
173  //delete theTrackExtra;
174  }
175 
177 
178  for (int k = 0; k < nTracks; k++) {
179 
180  const reco::TrackExtraRef theTrackExtraRef(ohTE,k);
181  (tracks->at(k)).setExtra(theTrackExtraRef);
182 
183  }
184 
185  e.put(std::move(tracks));
186 
187 }
188 
T getParameter(std::string const &) const
const unsigned int nTracks(const reco::Vertex &sv)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void setHits(TrackingRecHitRefProd const &prod, unsigned firstH, unsigned int nH)
PixelTracksProducer(const edm::ParameterSet &conf)
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
virtual void produce(edm::Event &ev, const edm::EventSetup &es)
TrackingRegionProducer * theRegionProducer
edm::EDGetTokenT< TrajectorySeedCollection > seedProducerToken
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
edm::EDGetTokenT< PixelTrackFilter > filterToken
std::pair< const_iterator, const_iterator > range
std::pair< reco::Track *, std::vector< const TrackingRecHit * > > TrackWithRecHits
int k[5][pyjets_maxn]
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:11
virtual std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &ev, const edm::EventSetup &es) const =0
const T & get() const
Definition: EventSetup.h:56
std::unique_ptr< reco::Track > run(const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region) const
Definition: PixelFitter.h:15
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135
bool appendHitPattern(const TrackingRecHit &hit, const TrackerTopology &ttopo)
append a single hit to the HitPattern
Definition: TrackBase.h:455
std::vector< TrackWithRecHits > TracksWithRecHits
edm::EDGetTokenT< PixelFitter > fitterToken
def move(src, dest)
Definition: eostools.py:510
T get(const Candidate &c)
Definition: component.h:55