CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
PixelTrackReconstruction Class Reference

#include <PixelTrackReconstruction.h>

Public Member Functions

void halt ()
 
void init (const edm::EventSetup &es)
 
 PixelTrackReconstruction (const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
 
void run (pixeltrackfitting::TracksWithTTRHs &tah, edm::Event &ev, const edm::EventSetup &es)
 
 ~PixelTrackReconstruction ()
 

Private Attributes

PixelTrackCleanertheCleaner
 
edm::ParameterSet theConfig
 
std::unique_ptr< PixelTrackFiltertheFilter
 
const PixelFittertheFitter
 
std::unique_ptr
< OrderedHitsGenerator
theGenerator
 
std::unique_ptr
< QuadrupletSeedMerger
theMerger_
 
std::unique_ptr
< TrackingRegionProducer
theRegionProducer
 
bool useClusterShape
 

Detailed Description

Definition at line 19 of file PixelTrackReconstruction.h.

Constructor & Destructor Documentation

PixelTrackReconstruction::PixelTrackReconstruction ( const edm::ParameterSet conf,
edm::ConsumesCollector &&  iC 
)

Definition at line 36 of file PixelTrackReconstruction.cc.

References SurfaceDeformationFactory::create(), edm::ParameterSet::exists(), reco_skim_cfg_mod::filterName, reco::get(), edm::ParameterSet::getParameter(), AlCaHLTBitMon_QueryRunRegistry::string, theConfig, theFilter, theGenerator, theMerger_, theRegionProducer, and useClusterShape.

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  useClusterShape = filterPSet.exists("useClusterShape");
60  }
61 
62  ParameterSet orderedPSet =
63  theConfig.getParameter<ParameterSet>("OrderedHitsFactoryPSet");
64  std::string orderedName = orderedPSet.getParameter<std::string>("ComponentName");
65  theGenerator.reset(OrderedHitsGeneratorFactory::get()->create( orderedName, orderedPSet, iC));
66 
67  ParameterSet regfactoryPSet = theConfig.getParameter<ParameterSet>("RegionFactoryPSet");
68  std::string regfactoryName = regfactoryPSet.getParameter<std::string>("ComponentName");
69  theRegionProducer.reset(TrackingRegionProducerFactory::get()->create(regfactoryName,regfactoryPSet, std::move(iC)));
70 }
T getParameter(std::string const &) const
tuple cfg
Definition: looper.py:237
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::unique_ptr< TrackingRegionProducer > theRegionProducer
std::unique_ptr< OrderedHitsGenerator > theGenerator
std::unique_ptr< PixelTrackFilter > theFilter
SurfaceDeformation * create(int type, const std::vector< double > &params)
std::unique_ptr< QuadrupletSeedMerger > theMerger_
T get(const Candidate &c)
Definition: component.h:55
PixelTrackReconstruction::~PixelTrackReconstruction ( )

Definition at line 72 of file PixelTrackReconstruction.cc.

References halt().

73 {
74  halt();
75 }

Member Function Documentation

void PixelTrackReconstruction::halt ( )

Definition at line 77 of file PixelTrackReconstruction.cc.

References theCleaner, and theFitter.

Referenced by PixelTrackProducer::endRun(), and ~PixelTrackReconstruction().

78 {
79  delete theFitter; theFitter=0;
80  delete theCleaner; theCleaner=0;
81 }
void PixelTrackReconstruction::init ( const edm::EventSetup es)

Definition at line 83 of file PixelTrackReconstruction.cc.

References reco::get(), edm::ParameterSet::getParameter(), AlCaHLTBitMon_QueryRunRegistry::string, theCleaner, theConfig, theFitter, and theMerger_.

Referenced by PixelTrackProducer::beginRun().

84 {
85 
86  ParameterSet fitterPSet = theConfig.getParameter<ParameterSet>("FitterPSet");
87  std::string fitterName = fitterPSet.getParameter<std::string>("ComponentName");
88  theFitter = PixelFitterFactory::get()->create( fitterName, fitterPSet);
89 
90  ParameterSet cleanerPSet = theConfig.getParameter<ParameterSet>("CleanerPSet");
91  std::string cleanerName = cleanerPSet.getParameter<std::string>("ComponentName");
92  if (cleanerName != "none") theCleaner = PixelTrackCleanerFactory::get()->create( cleanerName, cleanerPSet);
93 
94  if (theMerger_) {
95  theMerger_->update( es );
96  }
97 }
T getParameter(std::string const &) const
std::unique_ptr< QuadrupletSeedMerger > theMerger_
T get(const Candidate &c)
Definition: component.h:55
void PixelTrackReconstruction::run ( pixeltrackfitting::TracksWithTTRHs tah,
edm::Event ev,
const edm::EventSetup es 
)

Definition at line 99 of file PixelTrackReconstruction.cc.

References PixelTrackCleanerWrapper::clean(), edm::EventSetup::get(), edm::ESHandle< class >::product(), PixelFitter::run(), OrderedSeedingHits::size(), SeedingHitSet::size(), theCleaner, theFilter, theFitter, theGenerator, theMerger_, theRegionProducer, and useClusterShape.

Referenced by PixelTrackProducer::produce().

100 {
101  typedef std::vector<TrackingRegion* > Regions;
102  typedef Regions::const_iterator IR;
103  Regions regions = theRegionProducer->regions(ev,es);
104 
105  //Retrieve tracker topology from geometry
107  es.get<IdealGeometryRecord>().get(tTopoHand);
108  const TrackerTopology *tTopo=tTopoHand.product();
109 
110  if (theFilter) theFilter->update(ev, es);
111 
112  for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) {
113  const TrackingRegion & region = **ir;
114 
115  const OrderedSeedingHits & triplets = theGenerator->run(region,ev,es);
116  const OrderedSeedingHits &tuplets= (theMerger_==0)? triplets : theMerger_->mergeTriplets( triplets, es );
117 
118  unsigned int nTuplets = tuplets.size();
119  tracks.reserve(tracks.size()+nTuplets);
120  // producing tracks
121  for (unsigned int iTuplet = 0; iTuplet < nTuplets; ++iTuplet) {
122  const SeedingHitSet & tuplet = tuplets[iTuplet];
123 
124  std::vector<const TrackingRecHit *> hits;
125  for (unsigned int iHit = 0, nHits = tuplet.size(); iHit < nHits; ++iHit) {
126  hits.push_back( tuplet[iHit]->hit() );
127  }
128 
129  // fitting
130  reco::Track* track = theFitter->run( ev, es, hits, region);
131  if (!track) continue;
132 
133  if (theFilter) {
134  if ((useClusterShape && !(*theFilter)(track, hits, tTopo)) ||
135  (!useClusterShape && !(*theFilter)(track, hits))) {
136  delete track;
137  continue;
138  }
139  }
140 
141  // add tracks
142  tracks.push_back(TrackWithTTRHs(track, tuplet));
143  }
144  theGenerator->clear();
145  }
146 
147  // skip ovelrapped tracks
149 
150  // clean memory
151  for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) delete (*ir);
152 }
std::pair< reco::Track *, SeedingHitSet > TrackWithTTRHs
virtual reco::Track * run(const edm::EventSetup &es, const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region) const
Definition: PixelFitter.h:17
virtual unsigned int size() const =0
std::unique_ptr< TrackingRegionProducer > theRegionProducer
std::unique_ptr< OrderedHitsGenerator > theGenerator
std::unique_ptr< PixelTrackFilter > theFilter
tuple tracks
Definition: testEve_cfg.py:39
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:86
pixeltrackfitting::TracksWithTTRHs clean(const pixeltrackfitting::TracksWithTTRHs &initialT_TTRHs, const TrackerTopology *tTopo)
unsigned int size() const
Definition: SeedingHitSet.h:44
std::unique_ptr< QuadrupletSeedMerger > theMerger_

Member Data Documentation

PixelTrackCleaner* PixelTrackReconstruction::theCleaner
private

Definition at line 35 of file PixelTrackReconstruction.h.

Referenced by halt(), init(), and run().

edm::ParameterSet PixelTrackReconstruction::theConfig
private

Definition at line 32 of file PixelTrackReconstruction.h.

Referenced by init(), and PixelTrackReconstruction().

std::unique_ptr<PixelTrackFilter> PixelTrackReconstruction::theFilter
private

Definition at line 34 of file PixelTrackReconstruction.h.

Referenced by PixelTrackReconstruction(), and run().

const PixelFitter* PixelTrackReconstruction::theFitter
private

Definition at line 33 of file PixelTrackReconstruction.h.

Referenced by halt(), init(), and run().

std::unique_ptr<OrderedHitsGenerator> PixelTrackReconstruction::theGenerator
private

Definition at line 36 of file PixelTrackReconstruction.h.

Referenced by PixelTrackReconstruction(), and run().

std::unique_ptr<QuadrupletSeedMerger> PixelTrackReconstruction::theMerger_
private

Definition at line 38 of file PixelTrackReconstruction.h.

Referenced by init(), PixelTrackReconstruction(), and run().

std::unique_ptr<TrackingRegionProducer> PixelTrackReconstruction::theRegionProducer
private

Definition at line 37 of file PixelTrackReconstruction.h.

Referenced by PixelTrackReconstruction(), and run().

bool PixelTrackReconstruction::useClusterShape
private

Definition at line 39 of file PixelTrackReconstruction.h.

Referenced by PixelTrackReconstruction(), and run().