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)
 
void run (pixeltrackfitting::TracksWithTTRHs &tah, edm::Event &ev, const edm::EventSetup &es)
 
 ~PixelTrackReconstruction ()
 

Private Attributes

PixelTrackCleanertheCleaner
 
edm::ParameterSet theConfig
 
PixelTrackFiltertheFilter
 
const PixelFittertheFitter
 
OrderedHitsGeneratortheGenerator
 
QuadrupletSeedMergertheMerger_
 
TrackingRegionProducertheRegionProducer
 

Detailed Description

Definition at line 17 of file PixelTrackReconstruction.h.

Constructor & Destructor Documentation

PixelTrackReconstruction::PixelTrackReconstruction ( const edm::ParameterSet conf)

Definition at line 36 of file PixelTrackReconstruction.cc.

References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), AlCaHLTBitMon_QueryRunRegistry::string, theConfig, and theMerger_.

38 {
39  if ( cfg.exists("SeedMergerPSet") ) {
40  edm::ParameterSet mergerPSet = theConfig.getParameter<edm::ParameterSet>( "SeedMergerPSet" );
41  std::string seedmergerTTRHBuilderLabel = mergerPSet.getParameter<std::string>( "ttrhBuilderLabel" );
42  std::string seedmergerLayerListName = mergerPSet.getParameter<std::string>( "layerListName" );
43  bool seedmergerAddTriplets = mergerPSet.getParameter<bool>( "addRemainingTriplets" );
44  bool seedmergerMergeTriplets = mergerPSet.getParameter<bool>( "mergeTriplets" );
46  theMerger_->setMergeTriplets( seedmergerMergeTriplets );
47  theMerger_->setAddRemainingTriplets( seedmergerAddTriplets );
48  theMerger_->setTTRHBuilderLabel( seedmergerTTRHBuilderLabel );
49  theMerger_->setLayerListName( seedmergerLayerListName );
50  }
51 }
T getParameter(std::string const &) const
OrderedHitsGenerator * theGenerator
TrackingRegionProducer * theRegionProducer
QuadrupletSeedMerger * theMerger_
PixelTrackReconstruction::~PixelTrackReconstruction ( )

Definition at line 53 of file PixelTrackReconstruction.cc.

References halt().

54 {
55  halt();
56 }

Member Function Documentation

void PixelTrackReconstruction::halt ( )

Definition at line 58 of file PixelTrackReconstruction.cc.

References theCleaner, theFilter, theFitter, theGenerator, theMerger_, and theRegionProducer.

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

59 {
60  delete theFilter; theFilter=0;
61  delete theFitter; theFitter=0;
62  delete theCleaner; theCleaner=0;
63  delete theGenerator; theGenerator=0;
65  delete theMerger_; theMerger_=0;
66 }
OrderedHitsGenerator * theGenerator
TrackingRegionProducer * theRegionProducer
QuadrupletSeedMerger * theMerger_
void PixelTrackReconstruction::init ( const edm::EventSetup es)

Definition at line 68 of file PixelTrackReconstruction.cc.

References reco_skim_cfg_mod::filterName, reco::get(), edm::ParameterSet::getParameter(), AlCaHLTBitMon_QueryRunRegistry::string, theCleaner, theConfig, theFilter, theFitter, theGenerator, theMerger_, theRegionProducer, and QuadrupletSeedMerger::update().

Referenced by PixelTrackProducer::beginRun().

69 {
70  ParameterSet regfactoryPSet = theConfig.getParameter<ParameterSet>("RegionFactoryPSet");
71  std::string regfactoryName = regfactoryPSet.getParameter<std::string>("ComponentName");
72  theRegionProducer = TrackingRegionProducerFactory::get()->create(regfactoryName,regfactoryPSet);
73 
74  ParameterSet orderedPSet =
75  theConfig.getParameter<ParameterSet>("OrderedHitsFactoryPSet");
76  std::string orderedName = orderedPSet.getParameter<std::string>("ComponentName");
77  theGenerator = OrderedHitsGeneratorFactory::get()->create( orderedName, orderedPSet);
78 
79  ParameterSet fitterPSet = theConfig.getParameter<ParameterSet>("FitterPSet");
80  std::string fitterName = fitterPSet.getParameter<std::string>("ComponentName");
81  theFitter = PixelFitterFactory::get()->create( fitterName, fitterPSet);
82 
83  ParameterSet filterPSet = theConfig.getParameter<ParameterSet>("FilterPSet");
84  std::string filterName = filterPSet.getParameter<std::string>("ComponentName");
85  if (filterName != "none") {
86  theFilter = theConfig.getParameter<bool>("useFilterWithES") ?
87  PixelTrackFilterWithESFactory::get()->create( filterName, filterPSet, es) :
88  PixelTrackFilterFactory::get()->create( filterName, filterPSet);
89  }
90 
91  ParameterSet cleanerPSet = theConfig.getParameter<ParameterSet>("CleanerPSet");
92  std::string cleanerName = cleanerPSet.getParameter<std::string>("ComponentName");
93  if (cleanerName != "none") theCleaner = PixelTrackCleanerFactory::get()->create( cleanerName, cleanerPSet);
94 
95  if ( theMerger_ !=0 ) {
96  theMerger_->update( es );
97  }
98 }
T getParameter(std::string const &) const
OrderedHitsGenerator * theGenerator
TrackingRegionProducer * theRegionProducer
void update(const edm::EventSetup &)
T get(const Candidate &c)
Definition: component.h:56
QuadrupletSeedMerger * theMerger_
void PixelTrackReconstruction::run ( pixeltrackfitting::TracksWithTTRHs tah,
edm::Event ev,
const edm::EventSetup es 
)

Definition at line 100 of file PixelTrackReconstruction.cc.

References PixelTrackCleanerWrapper::clean(), OrderedHitsGenerator::clear(), edm::EventSetup::get(), QuadrupletSeedMerger::mergeTriplets(), edm::ESHandle< class >::product(), TrackingRegionProducer::regions(), OrderedHitsGenerator::run(), PixelFitter::run(), OrderedSeedingHits::size(), SeedingHitSet::size(), theCleaner, theFilter, theFitter, theGenerator, theMerger_, theRegionProducer, and PixelTrackFilter::update().

Referenced by PixelTrackProducer::produce().

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

Member Data Documentation

PixelTrackCleaner* PixelTrackReconstruction::theCleaner
private

Definition at line 32 of file PixelTrackReconstruction.h.

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

edm::ParameterSet PixelTrackReconstruction::theConfig
private

Definition at line 29 of file PixelTrackReconstruction.h.

Referenced by init(), and PixelTrackReconstruction().

PixelTrackFilter* PixelTrackReconstruction::theFilter
private

Definition at line 31 of file PixelTrackReconstruction.h.

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

const PixelFitter* PixelTrackReconstruction::theFitter
private

Definition at line 30 of file PixelTrackReconstruction.h.

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

OrderedHitsGenerator* PixelTrackReconstruction::theGenerator
private

Definition at line 33 of file PixelTrackReconstruction.h.

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

QuadrupletSeedMerger* PixelTrackReconstruction::theMerger_
private

Definition at line 35 of file PixelTrackReconstruction.h.

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

TrackingRegionProducer* PixelTrackReconstruction::theRegionProducer
private

Definition at line 34 of file PixelTrackReconstruction.h.

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