CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelTrackReconstruction.cc
Go to the documentation of this file.
2 
6 
9 
13 
16 
19 
23 
27 
29 
30 #include <vector>
31 
32 using namespace pixeltrackfitting;
33 using namespace ctfseeding;
34 using edm::ParameterSet;
35 
38  : theConfig(cfg), theFitter(0), theCleaner(0)
39 {
40  if ( cfg.exists("SeedMergerPSet") ) {
41  edm::ParameterSet mergerPSet = theConfig.getParameter<edm::ParameterSet>( "SeedMergerPSet" );
42  std::string seedmergerTTRHBuilderLabel = mergerPSet.getParameter<std::string>( "ttrhBuilderLabel" );
43  edm::ParameterSet seedmergerLayerList = mergerPSet.getParameter<edm::ParameterSet>( "layerList" );
44  bool seedmergerAddTriplets = mergerPSet.getParameter<bool>( "addRemainingTriplets" );
45  bool seedmergerMergeTriplets = mergerPSet.getParameter<bool>( "mergeTriplets" );
46  theMerger_.reset(new QuadrupletSeedMerger(seedmergerLayerList, iC));
47  theMerger_->setMergeTriplets( seedmergerMergeTriplets );
48  theMerger_->setAddRemainingTriplets( seedmergerAddTriplets );
49  theMerger_->setTTRHBuilderLabel( seedmergerTTRHBuilderLabel );
50  }
51 
52  ParameterSet filterPSet = theConfig.getParameter<ParameterSet>("FilterPSet");
53  std::string filterName = filterPSet.getParameter<std::string>("ComponentName");
54  if (filterName != "none") {
55  theFilter.reset(PixelTrackFilterFactory::get()->create( filterName, filterPSet, iC));
56  if(theConfig.exists("useFilterWithES")) {
57  edm::LogInfo("Obsolete") << "useFilterWithES parameter is obsolete and can be removed";
58  }
59  }
60 
61  ParameterSet orderedPSet =
62  theConfig.getParameter<ParameterSet>("OrderedHitsFactoryPSet");
63  std::string orderedName = orderedPSet.getParameter<std::string>("ComponentName");
64  theGenerator.reset(OrderedHitsGeneratorFactory::get()->create( orderedName, orderedPSet, iC));
65 
66  ParameterSet regfactoryPSet = theConfig.getParameter<ParameterSet>("RegionFactoryPSet");
67  std::string regfactoryName = regfactoryPSet.getParameter<std::string>("ComponentName");
68  theRegionProducer.reset(TrackingRegionProducerFactory::get()->create(regfactoryName,regfactoryPSet, std::move(iC)));
69 }
70 
72 {
73  halt();
74 }
75 
77 {
78  delete theFitter; theFitter=0;
79  delete theCleaner; theCleaner=0;
80 }
81 
83 {
84 
85  ParameterSet fitterPSet = theConfig.getParameter<ParameterSet>("FitterPSet");
86  std::string fitterName = fitterPSet.getParameter<std::string>("ComponentName");
87  theFitter = PixelFitterFactory::get()->create( fitterName, fitterPSet);
88 
89  ParameterSet cleanerPSet = theConfig.getParameter<ParameterSet>("CleanerPSet");
90  std::string cleanerName = cleanerPSet.getParameter<std::string>("ComponentName");
91  if (cleanerName != "none") theCleaner = PixelTrackCleanerFactory::get()->create( cleanerName, cleanerPSet);
92 
93  if (theMerger_) {
94  theMerger_->update( es );
95  }
96 }
97 
99 {
100  typedef std::vector<TrackingRegion* > Regions;
101  typedef Regions::const_iterator IR;
102  Regions regions = theRegionProducer->regions(ev,es);
103 
104  //Retrieve tracker topology from geometry
106  es.get<IdealGeometryRecord>().get(tTopoHand);
107  const TrackerTopology *tTopo=tTopoHand.product();
108 
109  if (theFilter) theFilter->update(ev, es);
110 
111  for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) {
112  const TrackingRegion & region = **ir;
113 
114  const OrderedSeedingHits & triplets = theGenerator->run(region,ev,es);
115  const OrderedSeedingHits &tuplets= (theMerger_==0)? triplets : theMerger_->mergeTriplets( triplets, es );
116 
117  unsigned int nTuplets = tuplets.size();
118  tracks.reserve(tracks.size()+nTuplets);
119  // producing tracks
120  for (unsigned int iTuplet = 0; iTuplet < nTuplets; ++iTuplet) {
121  const SeedingHitSet & tuplet = tuplets[iTuplet];
122 
123  std::vector<const TrackingRecHit *> hits;
124  for (unsigned int iHit = 0, nHits = tuplet.size(); iHit < nHits; ++iHit) {
125  hits.push_back( tuplet[iHit]->hit() );
126  }
127 
128  // fitting
129  reco::Track* track = theFitter->run( ev, es, hits, region);
130  if (!track) continue;
131 
132  // decide if track should be skipped according to filter
133  if (theFilter && !(*theFilter)(track, hits) ) {
134  delete track;
135  continue;
136  }
137 
138  // add tracks
139  tracks.push_back(TrackWithTTRHs(track, tuplet));
140  }
141  theGenerator->clear();
142  }
143 
144  // skip ovelrapped tracks
145  if (theCleaner) tracks = PixelTrackCleanerWrapper(theCleaner).clean(tracks,tTopo);
146 
147  // clean memory
148  for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) delete (*ir);
149 }
T getParameter(std::string const &) const
std::pair< reco::Track *, SeedingHitSet > TrackWithTTRHs
void init(const edm::EventSetup &es)
virtual reco::Track * run(const edm::EventSetup &es, const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region) const
Definition: PixelFitter.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
virtual unsigned int size() const =0
std::vector< TrackWithTTRHs > TracksWithTTRHs
std::unique_ptr< TrackingRegionProducer > theRegionProducer
std::unique_ptr< OrderedHitsGenerator > theGenerator
PixelTrackReconstruction(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
std::unique_ptr< PixelTrackFilter > theFilter
void run(pixeltrackfitting::TracksWithTTRHs &tah, edm::Event &ev, const edm::EventSetup &es)
tuple tracks
Definition: testEve_cfg.py:39
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
pixeltrackfitting::TracksWithTTRHs clean(const pixeltrackfitting::TracksWithTTRHs &initialT_TTRHs, const TrackerTopology *tTopo)
unsigned int size() const
Definition: SeedingHitSet.h:44
SurfaceDeformation * create(int type, const std::vector< double > &params)
std::unique_ptr< QuadrupletSeedMerger > theMerger_
T get(const Candidate &c)
Definition: component.h:55