CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
SeedFinderSelector Class Reference

#include <SeedFinderSelector.h>

Public Member Functions

void initEvent (const edm::Event &, const edm::EventSetup &)
 
SeedingLayerSetsBuilder::SeedingLayerId Layer_tuple (const FastTrackerRecHit *hit) const
 
bool pass (const std::vector< const FastTrackerRecHit *> &hits) const
 
 SeedFinderSelector (const edm::ParameterSet &, edm::ConsumesCollector &&)
 
void setTrackingRegion (const TrackingRegion *trackingRegion)
 
 ~SeedFinderSelector ()
 

Private Attributes

std::unique_ptr< CAHitQuadrupletGeneratorCAHitQuadGenerator_
 
std::unique_ptr< CAHitTripletGeneratorCAHitTriplGenerator_
 
const edm::EventSetupeventSetup_
 
const MagneticFieldfield_ = nullptr
 
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecordfieldESToken_
 
std::vector< unsigned > layerPairs_
 
const MeasurementTrackermeasurementTracker_
 
const edm::ESGetToken< MeasurementTracker, CkfComponentsRecordmeasurementTrackerESToken_
 
const std::string measurementTrackerLabel_
 
const MultipleScatteringParametrisationMakermsmaker_ = nullptr
 
const edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecordmsMakerESToken_
 
std::unique_ptr< MultiHitGeneratorFromPairAndLayersmultiHitGenerator_
 
std::unique_ptr< HitTripletGeneratorFromPairAndLayerspixelTripletGenerator_
 
std::unique_ptr< SeedingLayerSetsHitsseedingLayer
 
std::vector< SeedingLayerSetsBuilder::SeedingLayerIdseedingLayerIds
 
std::unique_ptr< SeedingLayerSetsBuilderseedingLayers_
 
const TrackerTopologytrackerTopology_ = nullptr
 
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcdtrackerTopologyESToken_
 
const TrackingRegiontrackingRegion_
 

Detailed Description

Definition at line 29 of file SeedFinderSelector.h.

Constructor & Destructor Documentation

◆ SeedFinderSelector()

SeedFinderSelector::SeedFinderSelector ( const edm::ParameterSet cfg,
edm::ConsumesCollector &&  consumesCollector 
)

Definition at line 25 of file SeedFinderSelector.cc.

References CAHitQuadGenerator_, CAHitTriplGenerator_, looper::cfg, Exception, get, edm::ParameterSet::getParameter(), HLT_2022v15_cff::InputTag, layerPairs_, multiHitGenerator_, pixelTripletGenerator_, seedingLayers_, and AlCaHLTBitMon_QueryRunRegistry::string.

26  : trackingRegion_(nullptr),
27  eventSetup_(nullptr),
28  measurementTracker_(nullptr),
29  measurementTrackerLabel_(cfg.getParameter<std::string>("measurementTracker")),
31  trackerTopologyESToken_(consumesCollector.esConsumes()),
32  fieldESToken_(consumesCollector.esConsumes()),
33  msMakerESToken_(consumesCollector.esConsumes()) {
34  if (cfg.exists("pixelTripletGeneratorFactory")) {
35  const edm::ParameterSet &tripletConfig = cfg.getParameter<edm::ParameterSet>("pixelTripletGeneratorFactory");
37  tripletConfig.getParameter<std::string>("ComponentName"), tripletConfig, consumesCollector);
38  }
39 
40  if (cfg.exists("MultiHitGeneratorFactory")) {
41  const edm::ParameterSet &tripletConfig = cfg.getParameter<edm::ParameterSet>("MultiHitGeneratorFactory");
43  tripletConfig.getParameter<std::string>("ComponentName"), tripletConfig, consumesCollector);
44  }
45 
46  if (cfg.exists("CAHitTripletGeneratorFactory")) {
47  const edm::ParameterSet &tripletConfig = cfg.getParameter<edm::ParameterSet>("CAHitTripletGeneratorFactory");
48  CAHitTriplGenerator_ = std::make_unique<CAHitTripletGenerator>(tripletConfig, consumesCollector);
49  seedingLayers_ = std::make_unique<SeedingLayerSetsBuilder>(
50  cfg,
51  consumesCollector,
52  //calling the new FastSim specific constructor to make SeedingLayerSetsHits pointer for triplet iterations
53  edm::InputTag("fastTrackerRecHits"));
54  layerPairs_ = cfg.getParameter<std::vector<unsigned>>("layerPairs"); //allowed layer pairs for CA triplets
55  }
56 
57  if (cfg.exists("CAHitQuadrupletGeneratorFactory")) {
58  const edm::ParameterSet &quadrupletConfig = cfg.getParameter<edm::ParameterSet>("CAHitQuadrupletGeneratorFactory");
59  CAHitQuadGenerator_ = std::make_unique<CAHitQuadrupletGenerator>(quadrupletConfig, consumesCollector);
60  //calling the new FastSim specific constructor to make SeedingLayerSetsHits pointer for quadruplet iterations
62  std::make_unique<SeedingLayerSetsBuilder>(cfg, consumesCollector, edm::InputTag("fastTrackerRecHits"));
63  layerPairs_ = cfg.getParameter<std::vector<unsigned>>("layerPairs"); //allowed layer pairs for CA quadruplets
64  }
65 
68  throw cms::Exception("FastSimTracking")
69  << "It is forbidden to specify together 'pixelTripletGeneratorFactory', 'CAHitTripletGeneratorFactory' and "
70  "'MultiHitGeneratorFactory' in configuration of SeedFinderSelection";
71  }
74  throw cms::Exception("FastSimTracking")
75  << "It is forbidden to specify 'CAHitQuadrupletGeneratorFactory' together with 'pixelTripletGeneratorFactory', "
76  "'CAHitTripletGeneratorFactory' or 'MultiHitGeneratorFactory' in configuration of SeedFinderSelection";
77  }
78 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::unique_ptr< CAHitQuadrupletGenerator > CAHitQuadGenerator_
const edm::EventSetup * eventSetup_
const edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > msMakerESToken_
std::unique_ptr< HitTripletGeneratorFromPairAndLayers > pixelTripletGenerator_
std::vector< unsigned > layerPairs_
const TrackingRegion * trackingRegion_
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > fieldESToken_
const edm::ESGetToken< MeasurementTracker, CkfComponentsRecord > measurementTrackerESToken_
std::unique_ptr< SeedingLayerSetsBuilder > seedingLayers_
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > trackerTopologyESToken_
const MeasurementTracker * measurementTracker_
#define get
const std::string measurementTrackerLabel_
std::unique_ptr< CAHitTripletGenerator > CAHitTriplGenerator_
std::unique_ptr< MultiHitGeneratorFromPairAndLayers > multiHitGenerator_

◆ ~SeedFinderSelector()

SeedFinderSelector::~SeedFinderSelector ( )

Definition at line 80 of file SeedFinderSelector.cc.

80 { ; }

Member Function Documentation

◆ initEvent()

void SeedFinderSelector::initEvent ( const edm::Event ev,
const edm::EventSetup es 
)

Definition at line 82 of file SeedFinderSelector.cc.

References CAHitQuadGenerator_, CAHitTriplGenerator_, makeMEIFBenchmarkPlots::ev, eventSetup_, field_, fieldESToken_, edm::EventSetup::getData(), measurementTracker_, measurementTrackerESToken_, msmaker_, msMakerESToken_, multiHitGenerator_, seedingLayer, seedingLayerIds, seedingLayers_, trackerTopology_, and trackerTopologyESToken_.

82  {
83  eventSetup_ = &es;
84 
89 
90  if (multiHitGenerator_) {
91  multiHitGenerator_->initES(es);
92  }
93 
94  //for CA triplet iterations
96  seedingLayer = seedingLayers_->makeSeedingLayerSetsHitsforFastSim(ev, es);
97  seedingLayerIds = seedingLayers_->layers();
98  CAHitTriplGenerator_->initEvent(ev, es);
99  }
100  //for CA quadruplet iterations
101  if (CAHitQuadGenerator_) {
102  seedingLayer = seedingLayers_->makeSeedingLayerSetsHitsforFastSim(ev, es);
103  seedingLayerIds = seedingLayers_->layers();
104  CAHitQuadGenerator_->initEvent(ev, es);
105  }
106 }
std::unique_ptr< CAHitQuadrupletGenerator > CAHitQuadGenerator_
const edm::EventSetup * eventSetup_
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
std::unique_ptr< SeedingLayerSetsHits > seedingLayer
const MagneticField * field_
const edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > msMakerESToken_
std::vector< SeedingLayerSetsBuilder::SeedingLayerId > seedingLayerIds
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > fieldESToken_
const edm::ESGetToken< MeasurementTracker, CkfComponentsRecord > measurementTrackerESToken_
std::unique_ptr< SeedingLayerSetsBuilder > seedingLayers_
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > trackerTopologyESToken_
const MeasurementTracker * measurementTracker_
const MultipleScatteringParametrisationMaker * msmaker_
std::unique_ptr< CAHitTripletGenerator > CAHitTriplGenerator_
const TrackerTopology * trackerTopology_
std::unique_ptr< MultiHitGeneratorFromPairAndLayers > multiHitGenerator_

◆ Layer_tuple()

SeedingLayerSetsBuilder::SeedingLayerId SeedFinderSelector::Layer_tuple ( const FastTrackerRecHit hit) const

Definition at line 326 of file SeedFinderSelector.cc.

References Barrel, GeomDetEnumerators::invalidDet, NegEndcap, GeomDetEnumerators::PixelBarrel, PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, GeomDetEnumerators::PixelEndcap, PosEndcap, TrackerTopology::pxbLayer(), TrackerTopology::pxfDisk(), TrackerTopology::pxfSide(), and trackerTopology_.

Referenced by pass().

326  {
329  int idLayer = 0;
330 
331  if ((hit->det()->geographicalId()).subdetId() == PixelSubdetector::PixelBarrel) {
333  side = TrackerDetSide::Barrel;
334  idLayer = trackerTopology_->pxbLayer(hit->det()->geographicalId());
335  } else if ((hit->det()->geographicalId()).subdetId() == PixelSubdetector::PixelEndcap) {
337  idLayer = trackerTopology_->pxfDisk(hit->det()->geographicalId());
338  if (trackerTopology_->pxfSide(hit->det()->geographicalId()) == 1) {
340  } else {
342  }
343  }
344  return std::make_tuple(subdet, side, idLayer);
345 }
unsigned int pxbLayer(const DetId &id) const
TrackerDetSide
Definition: TrackerDetSide.h:4
unsigned int pxfDisk(const DetId &id) const
unsigned int pxfSide(const DetId &id) const
const TrackerTopology * trackerTopology_

◆ pass()

bool SeedFinderSelector::pass ( const std::vector< const FastTrackerRecHit *> &  hits) const

Definition at line 108 of file SeedFinderSelector.cc.

References CAHitQuadGenerator_, CAHitTriplGenerator_, GeometricSearchTracker::detLayer(), HitPairGeneratorFromLayerPair::doublets(), eventSetup_, Exception, field_, trigObjTnPSource_cfi::filler, MeasurementTracker::geometricSearchTracker(), hfClusterShapes_cfi::hits, Layer_tuple(), layerPairs_, hgcalTBTopologyTester_cfi::layers, eostools::ls(), measurementTracker_, eostools::move(), msmaker_, multiHitGenerator_, pixelTracksMonitoring_cff::ntuplet, TrackingRegion::origin(), AlCaHLTBitMon_ParallelJobs::p, pixelTripletGenerator_, mps_fire::result, seedingLayer, seedingLayerIds, and trackingRegion_.

108  {
109  if (!measurementTracker_ || !eventSetup_) {
110  throw cms::Exception("FastSimTracking") << "ERROR: event not initialized";
111  }
112  if (!trackingRegion_) {
113  throw cms::Exception("FastSimTracking") << "ERROR: trackingRegion not set";
114  }
115 
116  // check the inner 2 hits
117  if (hits.size() < 2) {
118  throw cms::Exception("FastSimTracking") << "SeedFinderSelector::pass requires at least 2 hits";
119  }
120  const DetLayer *firstLayer =
121  measurementTracker_->geometricSearchTracker()->detLayer(hits[0]->det()->geographicalId());
122  const DetLayer *secondLayer =
123  measurementTracker_->geometricSearchTracker()->detLayer(hits[1]->det()->geographicalId());
124 
125  std::vector<BaseTrackerRecHit const *> firstHits{hits[0]};
126  std::vector<BaseTrackerRecHit const *> secondHits{hits[1]};
127 
128  const RecHitsSortedInPhi fhm(firstHits, trackingRegion_->origin(), firstLayer);
129  const RecHitsSortedInPhi shm(secondHits, trackingRegion_->origin(), secondLayer);
130 
131  HitDoublets result(fhm, shm);
133  *trackingRegion_, *firstLayer, *secondLayer, fhm, shm, *field_, *msmaker_, 0, result);
134 
135  if (result.empty()) {
136  return false;
137  }
138 
139  // check the inner 3 hits
141  if (hits.size() < 3) {
142  throw cms::Exception("FastSimTracking")
143  << "For the given configuration, SeedFinderSelector::pass requires at least 3 hits";
144  }
145  const DetLayer *thirdLayer =
146  measurementTracker_->geometricSearchTracker()->detLayer(hits[2]->det()->geographicalId());
147  std::vector<const DetLayer *> thirdLayerDetLayer(1, thirdLayer);
148  std::vector<BaseTrackerRecHit const *> thirdHits{hits[2]};
149  const RecHitsSortedInPhi thm(thirdHits, trackingRegion_->origin(), thirdLayer);
150  const RecHitsSortedInPhi *thmp = &thm;
151 
153  OrderedHitTriplets tripletresult;
154  pixelTripletGenerator_->hitTriplets(
155  *trackingRegion_, tripletresult, *eventSetup_, result, &thmp, thirdLayerDetLayer, 1);
156  return !tripletresult.empty();
157  } else if (multiHitGenerator_) {
158  OrderedMultiHits tripletresult;
159  multiHitGenerator_->hitTriplets(*trackingRegion_, tripletresult, result, &thmp, thirdLayerDetLayer, 1);
160  return !tripletresult.empty();
161  }
162  //new for Phase1
163  else if (CAHitTriplGenerator_) {
164  if (!seedingLayer)
165  throw cms::Exception("FastSimTracking") << "ERROR: SeedingLayers pointer not set for CATripletGenerator";
166 
168  //constructing IntermediateHitDoublets to be passed onto CAHitTripletGenerator::hitNtuplets()
170  const TrackingRegion &tr_ = *trackingRegion_;
171  auto filler = ihd.beginRegion(&tr_);
172 
173  //forming the SeedingLayerId of the hits
174  std::array<SeedingLayerSetsBuilder::SeedingLayerId, 3> hitPair;
175  hitPair[0] = Layer_tuple(hits[0]);
176  hitPair[1] = Layer_tuple(hits[1]);
177  hitPair[2] = Layer_tuple(hits[2]);
178 
179  //extracting the DetLayer of the hits
180  const DetLayer *fLayer =
181  measurementTracker_->geometricSearchTracker()->detLayer(hits[0]->det()->geographicalId());
182  const DetLayer *sLayer =
183  measurementTracker_->geometricSearchTracker()->detLayer(hits[1]->det()->geographicalId());
184  const DetLayer *tLayer =
185  measurementTracker_->geometricSearchTracker()->detLayer(hits[2]->det()->geographicalId());
186 
187  //converting FastTrackerRecHit hits to BaseTrackerRecHit
188  std::vector<BaseTrackerRecHit const *> fHits{hits[0]};
189  std::vector<BaseTrackerRecHit const *> sHits{hits[1]};
190  std::vector<BaseTrackerRecHit const *> tHits{hits[2]};
191 
192  //forming the SeedingLayerSet for the hit doublets
193  SeedingLayerSetsHits::SeedingLayerSet pairCandidate1, pairCandidate2;
196  for (const auto p : layerPairs_) {
197  pairCandidate = ls.slice(p, p + 2);
198  if (p == 0 && hitPair[0] == seedingLayerIds[pairCandidate[0].index()] &&
199  hitPair[1] == seedingLayerIds[pairCandidate[1].index()])
200  pairCandidate1 = pairCandidate;
201  if (p == 1 && hitPair[1] == seedingLayerIds[pairCandidate[0].index()] &&
202  hitPair[2] == seedingLayerIds[pairCandidate[1].index()])
203  pairCandidate2 = pairCandidate;
204  }
205  }
206 
207  //Important: hits of the layer to be added to LayerHitMapCache
208  auto &layerCache = filler.layerHitMapCache();
209 
210  //doublets for CA triplets from the allowed layer pair combinations:(0,1),(1,2) and storing in filler
211  const RecHitsSortedInPhi &firsthm = *layerCache.add(
212  pairCandidate1[0], std::make_unique<RecHitsSortedInPhi>(fHits, trackingRegion_->origin(), fLayer));
213  const RecHitsSortedInPhi &secondhm = *layerCache.add(
214  pairCandidate1[1], std::make_unique<RecHitsSortedInPhi>(sHits, trackingRegion_->origin(), sLayer));
215  HitDoublets res1(firsthm, secondhm);
217  *trackingRegion_, *fLayer, *sLayer, firsthm, secondhm, *field_, *msmaker_, 0, res1);
218  filler.addDoublets(pairCandidate1, std::move(res1));
219  const RecHitsSortedInPhi &thirdhm = *layerCache.add(
220  pairCandidate2[1], std::make_unique<RecHitsSortedInPhi>(tHits, trackingRegion_->origin(), tLayer));
221  HitDoublets res2(secondhm, thirdhm);
223  *trackingRegion_, *sLayer, *tLayer, secondhm, thirdhm, *field_, *msmaker_, 0, res2);
224  filler.addDoublets(pairCandidate2, std::move(res2));
225 
226  std::vector<OrderedHitSeeds> tripletresult;
227  tripletresult.resize(ihd.regionSize());
228  for (auto &ntuplet : tripletresult)
229  ntuplet.reserve(3);
230  //calling the function from the class, modifies tripletresult
231  CAHitTriplGenerator_->hitNtuplets(ihd, tripletresult, *seedingLayer);
232  return !tripletresult[0].empty();
233  }
234  }
235  //new for Phase1
236  if (CAHitQuadGenerator_) {
237  if (hits.size() < 4) {
238  throw cms::Exception("FastSimTracking")
239  << "For the given configuration, SeedFinderSelector::pass requires at least 4 hits";
240  }
241 
242  if (!seedingLayer)
243  throw cms::Exception("FastSimTracking") << "ERROR: SeedingLayers pointer not set for CAHitQuadrupletGenerator";
244 
246  //constructing IntermediateHitDoublets to be passed onto CAHitQuadrupletGenerator::hitNtuplets()
248  const TrackingRegion &tr_ = *trackingRegion_;
249  auto filler = ihd.beginRegion(&tr_);
250 
251  //forming the SeedingLayerId of the hits
252  std::array<SeedingLayerSetsBuilder::SeedingLayerId, 4> hitPair;
253  hitPair[0] = Layer_tuple(hits[0]);
254  hitPair[1] = Layer_tuple(hits[1]);
255  hitPair[2] = Layer_tuple(hits[2]);
256  hitPair[3] = Layer_tuple(hits[3]);
257 
258  //extracting the DetLayer of the hits
259  const DetLayer *fLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[0]->det()->geographicalId());
260  const DetLayer *sLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[1]->det()->geographicalId());
261  const DetLayer *tLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[2]->det()->geographicalId());
262  const DetLayer *frLayer = measurementTracker_->geometricSearchTracker()->detLayer(hits[3]->det()->geographicalId());
263 
264  //converting FastTrackerRecHit hits to BaseTrackerRecHit
265  std::vector<BaseTrackerRecHit const *> fHits{hits[0]};
266  std::vector<BaseTrackerRecHit const *> sHits{hits[1]};
267  std::vector<BaseTrackerRecHit const *> tHits{hits[2]};
268  std::vector<BaseTrackerRecHit const *> frHits{hits[3]};
269 
270  //forming the SeedingLayerSet for the hit doublets
271  SeedingLayerSetsHits::SeedingLayerSet pairCandidate1, pairCandidate2, pairCandidate3;
274  for (const auto p : layerPairs_) {
275  pairCandidate = ls.slice(p, p + 2);
276  if (p == 0 && hitPair[0] == seedingLayerIds[pairCandidate[0].index()] &&
277  hitPair[1] == seedingLayerIds[pairCandidate[1].index()])
278  pairCandidate1 = pairCandidate;
279  if (p == 1 && hitPair[1] == seedingLayerIds[pairCandidate[0].index()] &&
280  hitPair[2] == seedingLayerIds[pairCandidate[1].index()])
281  pairCandidate2 = pairCandidate;
282  if (p == 2 && hitPair[2] == seedingLayerIds[pairCandidate[0].index()] &&
283  hitPair[3] == seedingLayerIds[pairCandidate[1].index()])
284  pairCandidate3 = pairCandidate;
285  }
286  }
287 
288  //Important: hits of the layer to be added to LayerHitMapCache
289  auto &layerCache = filler.layerHitMapCache();
290 
291  //doublets for CA quadruplets from the allowed layer pair combinations:(0,1),(1,2),(2,3) and storing in filler
292  const RecHitsSortedInPhi &firsthm = *layerCache.add(
293  pairCandidate1[0], std::make_unique<RecHitsSortedInPhi>(fHits, trackingRegion_->origin(), fLayer));
294  const RecHitsSortedInPhi &secondhm = *layerCache.add(
295  pairCandidate1[1], std::make_unique<RecHitsSortedInPhi>(sHits, trackingRegion_->origin(), sLayer));
296  HitDoublets res1(firsthm, secondhm);
298  *trackingRegion_, *fLayer, *sLayer, firsthm, secondhm, *field_, *msmaker_, 0, res1);
299  filler.addDoublets(pairCandidate1, std::move(res1));
300  const RecHitsSortedInPhi &thirdhm = *layerCache.add(
301  pairCandidate2[1], std::make_unique<RecHitsSortedInPhi>(tHits, trackingRegion_->origin(), tLayer));
302  HitDoublets res2(secondhm, thirdhm);
304  *trackingRegion_, *sLayer, *tLayer, secondhm, thirdhm, *field_, *msmaker_, 0, res2);
305  filler.addDoublets(pairCandidate2, std::move(res2));
306  const RecHitsSortedInPhi &fourthhm = *layerCache.add(
307  pairCandidate3[1], std::make_unique<RecHitsSortedInPhi>(frHits, trackingRegion_->origin(), frLayer));
308  HitDoublets res3(thirdhm, fourthhm);
310  *trackingRegion_, *tLayer, *frLayer, thirdhm, fourthhm, *field_, *msmaker_, 0, res3);
311  filler.addDoublets(pairCandidate3, std::move(res3));
312 
313  std::vector<OrderedHitSeeds> quadrupletresult;
314  quadrupletresult.resize(ihd.regionSize());
315  for (auto &ntuplet : quadrupletresult)
316  ntuplet.reserve(4);
317  //calling the function from the class, modifies quadrupletresult
318  CAHitQuadGenerator_->hitNtuplets(ihd, quadrupletresult, *seedingLayer);
319  return !quadrupletresult[0].empty();
320  }
321 
322  return true;
323 }
std::unique_ptr< CAHitQuadrupletGenerator > CAHitQuadGenerator_
const edm::EventSetup * eventSetup_
std::unique_ptr< SeedingLayerSetsHits > seedingLayer
const MagneticField * field_
GlobalPoint const & origin() const
const DetLayer * detLayer(const DetId &id) const
obsolete method. Use idToLayer() instead.
SeedingLayerSetsBuilder::SeedingLayerId Layer_tuple(const FastTrackerRecHit *hit) const
std::unique_ptr< HitTripletGeneratorFromPairAndLayers > pixelTripletGenerator_
const GeometricSearchTracker * geometricSearchTracker() const
std::vector< SeedingLayerSetsBuilder::SeedingLayerId > seedingLayerIds
std::vector< unsigned > layerPairs_
const TrackingRegion * trackingRegion_
HitDoublets doublets(const TrackingRegion &reg, const edm::Event &ev, const edm::EventSetup &es, Layers layers)
def ls(path, rec=False)
Definition: eostools.py:349
const MeasurementTracker * measurementTracker_
const MultipleScatteringParametrisationMaker * msmaker_
std::unique_ptr< CAHitTripletGenerator > CAHitTriplGenerator_
def move(src, dest)
Definition: eostools.py:511
std::unique_ptr< MultiHitGeneratorFromPairAndLayers > multiHitGenerator_

◆ setTrackingRegion()

void SeedFinderSelector::setTrackingRegion ( const TrackingRegion trackingRegion)
inline

Definition at line 37 of file SeedFinderSelector.h.

References trackingRegion_.

37 { trackingRegion_ = trackingRegion; }
const TrackingRegion * trackingRegion_

Member Data Documentation

◆ CAHitQuadGenerator_

std::unique_ptr<CAHitQuadrupletGenerator> SeedFinderSelector::CAHitQuadGenerator_
private

Definition at line 58 of file SeedFinderSelector.h.

Referenced by initEvent(), pass(), and SeedFinderSelector().

◆ CAHitTriplGenerator_

std::unique_ptr<CAHitTripletGenerator> SeedFinderSelector::CAHitTriplGenerator_
private

Definition at line 57 of file SeedFinderSelector.h.

Referenced by initEvent(), pass(), and SeedFinderSelector().

◆ eventSetup_

const edm::EventSetup* SeedFinderSelector::eventSetup_
private

Definition at line 47 of file SeedFinderSelector.h.

Referenced by initEvent(), and pass().

◆ field_

const MagneticField* SeedFinderSelector::field_ = nullptr
private

Definition at line 55 of file SeedFinderSelector.h.

Referenced by initEvent(), and pass().

◆ fieldESToken_

const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> SeedFinderSelector::fieldESToken_
private

Definition at line 52 of file SeedFinderSelector.h.

Referenced by initEvent().

◆ layerPairs_

std::vector<unsigned> SeedFinderSelector::layerPairs_
private

Definition at line 61 of file SeedFinderSelector.h.

Referenced by pass(), and SeedFinderSelector().

◆ measurementTracker_

const MeasurementTracker* SeedFinderSelector::measurementTracker_
private

Definition at line 48 of file SeedFinderSelector.h.

Referenced by initEvent(), and pass().

◆ measurementTrackerESToken_

const edm::ESGetToken<MeasurementTracker, CkfComponentsRecord> SeedFinderSelector::measurementTrackerESToken_
private

Definition at line 50 of file SeedFinderSelector.h.

Referenced by initEvent().

◆ measurementTrackerLabel_

const std::string SeedFinderSelector::measurementTrackerLabel_
private

Definition at line 49 of file SeedFinderSelector.h.

◆ msmaker_

const MultipleScatteringParametrisationMaker* SeedFinderSelector::msmaker_ = nullptr
private

Definition at line 56 of file SeedFinderSelector.h.

Referenced by initEvent(), and pass().

◆ msMakerESToken_

const edm::ESGetToken<MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord> SeedFinderSelector::msMakerESToken_
private

Definition at line 53 of file SeedFinderSelector.h.

Referenced by initEvent().

◆ multiHitGenerator_

std::unique_ptr<MultiHitGeneratorFromPairAndLayers> SeedFinderSelector::multiHitGenerator_
private

Definition at line 45 of file SeedFinderSelector.h.

Referenced by initEvent(), pass(), and SeedFinderSelector().

◆ pixelTripletGenerator_

std::unique_ptr<HitTripletGeneratorFromPairAndLayers> SeedFinderSelector::pixelTripletGenerator_
private

Definition at line 44 of file SeedFinderSelector.h.

Referenced by pass(), and SeedFinderSelector().

◆ seedingLayer

std::unique_ptr<SeedingLayerSetsHits> SeedFinderSelector::seedingLayer
private

Definition at line 60 of file SeedFinderSelector.h.

Referenced by initEvent(), and pass().

◆ seedingLayerIds

std::vector<SeedingLayerSetsBuilder::SeedingLayerId> SeedFinderSelector::seedingLayerIds
private

Definition at line 62 of file SeedFinderSelector.h.

Referenced by initEvent(), and pass().

◆ seedingLayers_

std::unique_ptr<SeedingLayerSetsBuilder> SeedFinderSelector::seedingLayers_
private

Definition at line 59 of file SeedFinderSelector.h.

Referenced by initEvent(), and SeedFinderSelector().

◆ trackerTopology_

const TrackerTopology* SeedFinderSelector::trackerTopology_ = nullptr
private

Definition at line 54 of file SeedFinderSelector.h.

Referenced by initEvent(), and Layer_tuple().

◆ trackerTopologyESToken_

const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> SeedFinderSelector::trackerTopologyESToken_
private

Definition at line 51 of file SeedFinderSelector.h.

Referenced by initEvent().

◆ trackingRegion_

const TrackingRegion* SeedFinderSelector::trackingRegion_
private

Definition at line 46 of file SeedFinderSelector.h.

Referenced by pass(), and setTrackingRegion().