test
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 37 of file PixelTrackReconstruction.cc.

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

40 {
41  if ( cfg.exists("SeedMergerPSet") ) {
42  edm::ParameterSet mergerPSet = theConfig.getParameter<edm::ParameterSet>( "SeedMergerPSet" );
43  std::string seedmergerTTRHBuilderLabel = mergerPSet.getParameter<std::string>( "ttrhBuilderLabel" );
44  edm::ParameterSet seedmergerLayerList = mergerPSet.getParameter<edm::ParameterSet>( "layerList" );
45  bool seedmergerAddTriplets = mergerPSet.getParameter<bool>( "addRemainingTriplets" );
46  bool seedmergerMergeTriplets = mergerPSet.getParameter<bool>( "mergeTriplets" );
47  theMerger_.reset(new QuadrupletSeedMerger(seedmergerLayerList, iC));
48  theMerger_->setMergeTriplets( seedmergerMergeTriplets );
49  theMerger_->setAddRemainingTriplets( seedmergerAddTriplets );
50  theMerger_->setTTRHBuilderLabel( seedmergerTTRHBuilderLabel );
51  }
52 
53  ParameterSet filterPSet = theConfig.getParameter<ParameterSet>("FilterPSet");
54  std::string filterName = filterPSet.getParameter<std::string>("ComponentName");
55  if (filterName != "none") {
56  theFilter.reset(PixelTrackFilterFactory::get()->create( filterName, filterPSet, iC));
57  if(theConfig.exists("useFilterWithES")) {
58  edm::LogInfo("Obsolete") << "useFilterWithES parameter is obsolete and can be removed";
59  }
60  useClusterShape = filterPSet.exists("useClusterShape");
61  }
62 
63  ParameterSet orderedPSet =
64  theConfig.getParameter<ParameterSet>("OrderedHitsFactoryPSet");
65  std::string orderedName = orderedPSet.getParameter<std::string>("ComponentName");
66  theGenerator.reset(OrderedHitsGeneratorFactory::get()->create( orderedName, orderedPSet, iC));
67 
68  ParameterSet regfactoryPSet = theConfig.getParameter<ParameterSet>("RegionFactoryPSet");
69  std::string regfactoryName = regfactoryPSet.getParameter<std::string>("ComponentName");
70  theRegionProducer.reset(TrackingRegionProducerFactory::get()->create(regfactoryName,regfactoryPSet, std::move(iC)));
71 }
T getParameter(std::string const &) const
tuple cfg
Definition: looper.py:293
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::unique_ptr< TrackingRegionProducer > theRegionProducer
std::unique_ptr< OrderedHitsGenerator > theGenerator
def move
Definition: eostools.py:510
std::unique_ptr< PixelTrackFilter > theFilter
std::unique_ptr< QuadrupletSeedMerger > theMerger_
T get(const Candidate &c)
Definition: component.h:55
PixelTrackReconstruction::~PixelTrackReconstruction ( )

Definition at line 73 of file PixelTrackReconstruction.cc.

References halt().

74 {
75  halt();
76 }

Member Function Documentation

void PixelTrackReconstruction::halt ( )

Definition at line 78 of file PixelTrackReconstruction.cc.

References theCleaner, and theFitter.

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

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

Definition at line 84 of file PixelTrackReconstruction.cc.

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

Referenced by PixelTrackProducer::beginRun().

85 {
86 
87  ParameterSet fitterPSet = theConfig.getParameter<ParameterSet>("FitterPSet");
88  std::string fitterName = fitterPSet.getParameter<std::string>("ComponentName");
89  theFitter = PixelFitterFactory::get()->create( fitterName, fitterPSet);
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_) {
96  theMerger_->update( es );
97  }
98 }
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 
)

FIXME at some point we need to migrate the fitter...

Definition at line 100 of file PixelTrackReconstruction.cc.

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

Referenced by PixelTrackProducer::produce().

101 {
102  typedef std::vector<std::unique_ptr<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<TrackerTopologyRcd>().get(tTopoHand);
109  const TrackerTopology *tTopo=tTopoHand.product();
110 
111  if (theFilter) theFilter->update(ev, es);
112 
113  std::vector<const TrackingRecHit *> hits;hits.reserve(4);
114  for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) {
115  const TrackingRegion & region = **ir;
116 
117  const OrderedSeedingHits & triplets = theGenerator->run(region,ev,es);
118  const OrderedSeedingHits &tuplets= (theMerger_==0)? triplets : theMerger_->mergeTriplets( triplets, es );
119 
120  unsigned int nTuplets = tuplets.size();
121  tracks.reserve(tracks.size()+nTuplets);
122  // producing tracks
123  for (unsigned int iTuplet = 0; iTuplet < nTuplets; ++iTuplet) {
124  const SeedingHitSet & tuplet = tuplets[iTuplet];
125 
127  auto nHits = tuplet.size(); hits.resize(nHits);
128  for (unsigned int iHit = 0; iHit < nHits; ++iHit) hits[iHit] = tuplet[iHit];
129 
130  // fitting
131  reco::Track* track = theFitter->run( ev, es, hits, region);
132  if (!track) continue;
133 
134  if (theFilter) {
135  if ((useClusterShape && !(*theFilter)(track, hits, tTopo)) ||
136  (!useClusterShape && !(*theFilter)(track, hits))) {
137  delete track;
138  continue;
139  }
140  }
141 
142  // add tracks
143  tracks.emplace_back(track, tuplet);
144  }
145  theGenerator->clear();
146  }
147 
148  // skip ovelrapped tracks
149 
150  if (theCleaner) {
151  if (theCleaner->fast)
152  theCleaner->cleanTracks(tracks,tTopo);
153  else
155  }
156 }
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
virtual TracksWithRecHits cleanTracks(const TracksWithRecHits &tracksWithRecHits, const TrackerTopology *tTopo)
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:56
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:46
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().