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 Types | Private Member Functions | Private Attributes
PixelTripletLargeTipGenerator Class Reference

#include <PixelTripletLargeTipGenerator.h>

Inheritance diagram for PixelTripletLargeTipGenerator:
HitTripletGeneratorFromPairAndLayers HitTripletGenerator OrderedHitsGenerator

Public Member Functions

virtual void hitTriplets (const TrackingRegion &region, OrderedHitTriplets &trs, const edm::Event &ev, const edm::EventSetup &es)
 
void init (const HitPairGenerator &pairs, LayerCacheType *layerCache) override
 
const HitPairGeneratorpairGenerator () const
 
 PixelTripletLargeTipGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
 
void setSeedingLayers (SeedingLayerSetsHits::SeedingLayerSet pairLayers, std::vector< SeedingLayerSetsHits::SeedingLayer > thirdLayers) override
 
virtual ~PixelTripletLargeTipGenerator ()
 
- Public Member Functions inherited from HitTripletGeneratorFromPairAndLayers
virtual ~HitTripletGeneratorFromPairAndLayers ()
 
- Public Member Functions inherited from HitTripletGenerator
virtual void clear () final
 
 HitTripletGenerator (unsigned int size=500)
 
 HitTripletGenerator (HitTripletGenerator const &other)
 
virtual void hitTriplets (const TrackingRegion &reg, OrderedHitTriplets &prs, const edm::EventSetup &es)
 
virtual const OrderedHitTripletsrun (const TrackingRegion &region, const edm::Event &ev, const edm::EventSetup &es) final
 
virtual ~HitTripletGenerator ()
 
- Public Member Functions inherited from OrderedHitsGenerator
 OrderedHitsGenerator ()
 
virtual ~OrderedHitsGenerator ()
 

Private Types

typedef
CombinedHitTripletGenerator::LayerCacheType 
LayerCacheType
 

Private Member Functions

bool checkPhiInRange (float phi, float phi1, float phi2) const
 
std::pair< float, float > mergePhiRanges (const std::pair< float, float > &r1, const std::pair< float, float > &r2) const
 

Private Attributes

float dphi
 
float extraHitRPhitolerance
 
float extraHitRZtolerance
 
LayerCacheTypetheLayerCache
 
std::vector
< SeedingLayerSetsHits::SeedingLayer
theLayers
 
HitPairGeneratorthePairGenerator
 
bool useBend
 
bool useFixedPreFiltering
 
bool useMScat
 

Additional Inherited Members

- Public Types inherited from HitTripletGeneratorFromPairAndLayers
typedef LayerHitMapCache LayerCacheType
 
- Public Attributes inherited from OrderedHitsGenerator
unsigned int theMaxElement
 

Detailed Description

A HitTripletGenerator from HitPairGenerator and vector of Layers. The HitPairGenerator provides a set of hit pairs. For each pair the search for compatible hit(s) is done among provided Layers

Definition at line 21 of file PixelTripletLargeTipGenerator.h.

Member Typedef Documentation

Definition at line 23 of file PixelTripletLargeTipGenerator.h.

Constructor & Destructor Documentation

PixelTripletLargeTipGenerator::PixelTripletLargeTipGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)

Definition at line 43 of file PixelTripletLargeTipGenerator.cc.

References dphi, edm::ParameterSet::getParameter(), OrderedHitsGenerator::theMaxElement, and useFixedPreFiltering.

44  : thePairGenerator(0),
45  theLayerCache(0),
46  useFixedPreFiltering(cfg.getParameter<bool>("useFixedPreFiltering")),
47  extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
48  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
49  useMScat(cfg.getParameter<bool>("useMultScattering")),
50  useBend(cfg.getParameter<bool>("useBending"))
51 { theMaxElement=cfg.getParameter<unsigned int>("maxElement");
53  dphi = cfg.getParameter<double>("phiPreFiltering");
54 }
T getParameter(std::string const &) const
virtual PixelTripletLargeTipGenerator::~PixelTripletLargeTipGenerator ( )
inlinevirtual

Definition at line 28 of file PixelTripletLargeTipGenerator.h.

References thePairGenerator.

28 { delete thePairGenerator; }

Member Function Documentation

bool PixelTripletLargeTipGenerator::checkPhiInRange ( float  phi,
float  phi1,
float  phi2 
) const
private

Definition at line 364 of file PixelTripletLargeTipGenerator.cc.

References M_PI.

Referenced by hitTriplets().

365 { while (phi > phi2) phi -= 2. * M_PI;
366  while (phi < phi1) phi += 2. * M_PI;
367  return phi <= phi2;
368 }
#define M_PI
Definition: DDAxes.h:10
void PixelTripletLargeTipGenerator::hitTriplets ( const TrackingRegion region,
OrderedHitTriplets trs,
const edm::Event ev,
const edm::EventSetup es 
)
virtual

Implements HitTripletGenerator.

Definition at line 83 of file PixelTripletLargeTipGenerator.cc.

References angle(), KDTreeLinkerAlgo< DATA >::build(), checkPhiInRange(), PixelRecoUtilities::curvature(), HitDoublets::detLayer(), HitPairGenerator::doublets(), dphi, relativeConstraints::empty, ev, extraHitRPhitolerance, extraHitRZtolerance, f, plotBeamSpotDB::first, fnSigmaRZ, Geom::fpi(), Geom::ftwoPi(), edm::EventSetup::get(), i, HitDoublets::inner, DetLayer::isBarrel(), geometryCSVtoXML::line, bookConverter::max, ThirdHitPredictionFromCircle::HelixRZ::maxCurvature(), mergePhiRanges(), min(), nSigmaPhi, nSigmaRZ, TrackingRegion::originRBound(), HitDoublets::outer, p3, PV3DBase< T, PVType, FrameType >::perp(), edm::ESHandle< class >::product(), TrackingRegion::ptMin(), Basic2DVector< T >::r(), CosmicsPD_Skims::radius, HLT_25ns14e33_v3_cff::region, KDTreeLinkerAlgo< DATA >::search(), DetLayer::seqNum(), OrderedHitTriplets::size(), findQualityFiles::size, swap(), theLayers, OrderedHitsGenerator::theMaxElement, thePairGenerator, patCandidatesForDimuonsSequences_cff::tracker, useBend, useFixedPreFiltering, useMScat, findQualityFiles::v, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

87 {
89  es.get<TrackerDigiGeometryRecord>().get(tracker);
90 
91  //Retrieve tracker topology from geometry
93  es.get<IdealGeometryRecord>().get(tTopoHand);
94  const TrackerTopology *tTopo=tTopoHand.product();
95 
96  auto const & doublets = thePairGenerator->doublets(region,ev,es);
97 
98  if (doublets.empty()) return;
99 
100  auto outSeq = doublets.detLayer(HitDoublets::outer)->seqNum();
101 
102 
103  int size = theLayers.size();
104 
105 
106  using NodeInfo = KDTreeNodeInfo<unsigned int>;
107  std::vector<NodeInfo > layerTree; // re-used throughout
108  std::vector<unsigned int> foundNodes; // re-used throughout
109  foundNodes.reserve(100);
110  #ifdef __clang__
111  std::vector<KDTreeLinkerAlgo<unsigned int>> hitTree(size);
112  std::vector<LayerRZPredictions> mapPred(size);
113  #else
115  LayerRZPredictions mapPred[size];
116  #endif
117 
118  float rzError[size]; //save maximum errors
119  float maxphi = Geom::ftwoPi(), minphi = -maxphi; //increase to cater for any range
120 
121 
122  const RecHitsSortedInPhi * thirdHitMap[size];
123 
124  for(int il = 0; il < size; il++) {
125  thirdHitMap[il] = &(*theLayerCache)(theLayers[il], region, ev, es);
126  auto const & hits = *thirdHitMap[il];
127 
128  const DetLayer *layer = theLayers[il].detLayer();
129  LayerRZPredictions &predRZ = mapPred[il];
130  predRZ.line.initLayer(layer);
131  predRZ.helix1.initLayer(layer);
132  predRZ.helix2.initLayer(layer);
133  predRZ.line.initTolerance(extraHitRZtolerance);
134  predRZ.helix1.initTolerance(extraHitRZtolerance);
135  predRZ.helix2.initTolerance(extraHitRZtolerance);
136  predRZ.rzPositionFixup = MatchedHitRZCorrectionFromBending(layer,tTopo);
137 
138  layerTree.clear();
139  float minv=999999.0; float maxv = -999999.0; // Initialise to extreme values in case no hits
140  float maxErr=0.0f;
141  for (unsigned int i=0; i!=hits.size(); ++i) {
142  auto angle = hits.phi(i);
143  auto v = hits.gv(i);
144  //use (phi,r) for endcaps rather than (phi,z)
145  minv = std::min(minv,v); maxv = std::max(maxv,v);
146  float myerr = hits.dv[i];
147  maxErr = std::max(maxErr,myerr);
148  layerTree.emplace_back(i, angle, v); // save it
149  if (angle < 0) // wrap all points in phi
150  { layerTree.emplace_back(i, angle+Geom::ftwoPi(), v);}
151  else
152  { layerTree.emplace_back(i, angle-Geom::ftwoPi(), v);}
153  }
154  KDTreeBox phiZ(minphi, maxphi, minv-0.01f, maxv+0.01f); // declare our bounds
155  //add fudge factors in case only one hit and also for floating-point inaccuracy
156  hitTree[il].build(layerTree, phiZ); // make KDtree
157  rzError[il] = maxErr; //save error
158  }
159 
160  double curv = PixelRecoUtilities::curvature(1. / region.ptMin(), es);
161 
162  for (std::size_t ip =0; ip!=doublets.size(); ip++) {
163  auto xi = doublets.x(ip,HitDoublets::inner);
164  auto yi = doublets.y(ip,HitDoublets::inner);
165  auto zi = doublets.z(ip,HitDoublets::inner);
166  // auto rvi = doublets.rv(ip,HitDoublets::inner);
167  auto xo = doublets.x(ip,HitDoublets::outer);
168  auto yo = doublets.y(ip,HitDoublets::outer);
169  auto zo = doublets.z(ip,HitDoublets::outer);
170  // auto rvo = doublets.rv(ip,HitDoublets::outer);
171  GlobalPoint gp1(xi,yi,zi);
172  GlobalPoint gp2(xo,yo,zo);
173 
174  PixelRecoLineRZ line(gp1, gp2);
175  PixelRecoPointRZ point2(gp2.perp(), zo);
176  ThirdHitPredictionFromCircle predictionRPhi(gp1, gp2, extraHitRPhitolerance);
177 
178  Range generalCurvature = predictionRPhi.curvature(region.originRBound());
179  if (!intersect(generalCurvature, Range(-curv, curv))) continue;
180 
181  for(int il = 0; il < size; il++) {
182  if (hitTree[il].empty()) continue; // Don't bother if no hits
183  const DetLayer *layer = theLayers[il].detLayer();
184  bool barrelLayer = layer->isBarrel();
185 
186  Range curvature = generalCurvature;
187  ThirdHitCorrection correction(es, region.ptMin(), layer, line, point2, outSeq, useMScat);
188 
189  LayerRZPredictions &predRZ = mapPred[il];
190  predRZ.line.initPropagator(&line);
191 
192  Range rzRange;
193  if (useBend) {
194  // For the barrel region:
195  // swiping the helix passing through the two points across from
196  // negative to positive bending, can give us a sort of U-shaped
197  // projection onto the phi-z (barrel) or r-z plane (forward)
198  // so we checking minimum/maximum of all three possible extrema
199  //
200  // For the endcap region:
201  // Checking minimum/maximum radius of the helix projection
202  // onto an endcap plane, here we have to guard against
203  // looping tracks, when phi(delta z) gets out of control.
204  // HelixRZ::rAtZ should not follow looping tracks, but clamp
205  // to the minimum reachable r with the next-best lower |curvature|.
206  // So same procedure as for the barrel region can be applied.
207  //
208  // In order to avoid looking for potential looping tracks at all
209  // we also clamp the allowed curvature range for this layer,
210  // and potentially fail the layer entirely
211 
212  if (!barrelLayer) {
213  Range z3s = predRZ.line.detRange();
214  double z3 = z3s.first < 0 ? std::max(z3s.first, z3s.second)
215  : std::min(z3s.first, z3s.second);
216  double maxCurvature = HelixRZ::maxCurvature(&predictionRPhi,
217  gp1.z(), gp2.z(), z3);
218  if (!intersect(curvature, Range(-maxCurvature, maxCurvature)))
219  continue;
220  }
221 
222  HelixRZ helix1(&predictionRPhi, gp1.z(), gp2.z(), curvature.first);
223  HelixRZ helix2(&predictionRPhi, gp1.z(), gp2.z(), curvature.second);
224 
225  predRZ.helix1.initPropagator(&helix1);
226  predRZ.helix2.initPropagator(&helix2);
227 
228  Range rzRanges[2] = { predRZ.helix1(), predRZ.helix2() };
229  predRZ.helix1.initPropagator(nullptr);
230  predRZ.helix2.initPropagator(nullptr);
231 
232  rzRange.first = std::min(rzRanges[0].first, rzRanges[1].first);
233  rzRange.second = std::max(rzRanges[0].second, rzRanges[1].second);
234 
235  // if the allowed curvatures include a straight line,
236  // this can give us another extremum for allowed r/z
237  if (curvature.first * curvature.second < 0.0) {
238  Range rzLineRange = predRZ.line();
239  rzRange.first = std::min(rzRange.first, rzLineRange.first);
240  rzRange.second = std::max(rzRange.second, rzLineRange.second);
241  }
242  } else {
243  rzRange = predRZ.line();
244  }
245 
246  if (rzRange.first >= rzRange.second)
247  continue;
248 
249  correction.correctRZRange(rzRange);
250 
251  Range phiRange;
252  if (useFixedPreFiltering) {
253  float phi0 = doublets.phi(ip,HitDoublets::outer);
254  phiRange = Range(phi0 - dphi, phi0 + dphi);
255  } else {
256  Range radius;
257 
258  if (barrelLayer) {
259  radius = predRZ.line.detRange();
260  if (!intersect(rzRange, predRZ.line.detSize()))
261  continue;
262  } else {
263  radius = rzRange;
264  if (!intersect(radius, predRZ.line.detSize()))
265  continue;
266  }
267 
268  Range rPhi1 = predictionRPhi(curvature, radius.first);
269  Range rPhi2 = predictionRPhi(curvature, radius.second);
270  correction.correctRPhiRange(rPhi1);
271  correction.correctRPhiRange(rPhi2);
272  rPhi1.first /= radius.first;
273  rPhi1.second /= radius.first;
274  rPhi2.first /= radius.second;
275  rPhi2.second /= radius.second;
276  phiRange = mergePhiRanges(rPhi1, rPhi2);
277  }
278 
279  foundNodes.clear(); // Now recover hits in bounding box...
280  float prmin=phiRange.min(), prmax=phiRange.max(); //get contiguous range
281  if ((prmax-prmin) > Geom::ftwoPi())
282  { prmax=Geom::fpi(); prmin = -Geom::fpi();}
283  else
284  { while (prmax>maxphi) { prmin -= Geom::ftwoPi(); prmax -= Geom::ftwoPi();}
285  while (prmin<minphi) { prmin += Geom::ftwoPi(); prmax += Geom::ftwoPi();}
286  // This needs range -twoPi to +twoPi to work
287  }
288  if (barrelLayer) {
289  Range regMax = predRZ.line.detRange();
290  Range regMin = predRZ.line(regMax.min());
291  regMax = predRZ.line(regMax.max());
292  correction.correctRZRange(regMin);
293  correction.correctRZRange(regMax);
294  if (regMax.min() < regMin.min()) { swap(regMax, regMin);}
295  KDTreeBox phiZ(prmin, prmax,
296  regMin.min()-fnSigmaRZ*rzError[il],
297  regMax.max()+fnSigmaRZ*rzError[il]);
298  hitTree[il].search(phiZ, foundNodes);
299  }
300  else {
301  KDTreeBox phiZ(prmin, prmax,
302  rzRange.min()-fnSigmaRZ*rzError[il],
303  rzRange.max()+fnSigmaRZ*rzError[il]);
304  hitTree[il].search(phiZ, foundNodes);
305  }
306 
307  MatchedHitRZCorrectionFromBending l2rzFixup(doublets.hit(ip,HitDoublets::outer)->det()->geographicalId(), tTopo);
308  MatchedHitRZCorrectionFromBending l3rzFixup = predRZ.rzPositionFixup;
309 
310  thirdHitMap[il] = &(*theLayerCache)(theLayers[il], region, ev, es);
311  auto const & hits = *thirdHitMap[il];
312  for (auto KDdata : foundNodes) {
313  GlobalPoint p3 = hits.gp(KDdata);
314  double p3_r = p3.perp();
315  double p3_z = p3.z();
316  float p3_phi = hits.phi(KDdata);
317 
318  Range rangeRPhi = predictionRPhi(curvature, p3_r);
319  correction.correctRPhiRange(rangeRPhi);
320 
321  float ir = 1.f/p3_r;
322  float phiErr = nSigmaPhi * hits.drphi[KDdata]*ir;
323  if (!checkPhiInRange(p3_phi, rangeRPhi.first*ir-phiErr, rangeRPhi.second*ir+phiErr))
324  continue;
325 
326  Basic2DVector<double> thc(p3.x(), p3.y());
327 
328  auto curv_ = predictionRPhi.curvature(thc);
329  double p2_r = point2.r(); double p2_z = point2.z(); // they will be modified!
330 
331  l2rzFixup(predictionRPhi, curv_, *doublets.hit(ip,HitDoublets::outer), p2_r, p2_z, tTopo);
332  l3rzFixup(predictionRPhi, curv_, *hits.theHits[KDdata].hit(), p3_r, p3_z, tTopo);
333 
334  Range rangeRZ;
335  if (useBend) {
336  HelixRZ updatedHelix(&predictionRPhi, gp1.z(), p2_z, curv_);
337  rangeRZ = predRZ.helix1(barrelLayer ? p3_r : p3_z, updatedHelix);
338  } else {
339  float tIP = predictionRPhi.transverseIP(thc);
340  PixelRecoPointRZ updatedPoint2(p2_r, p2_z);
341  PixelRecoLineRZ updatedLine(line.origin(), point2, tIP);
342  rangeRZ = predRZ.line(barrelLayer ? p3_r : p3_z, line);
343  }
344  correction.correctRZRange(rangeRZ);
345 
346  double err = nSigmaRZ * hits.dv[KDdata];
347 
348  rangeRZ.first -= err, rangeRZ.second += err;
349 
350  if (!rangeRZ.inside(barrelLayer ? p3_z : p3_r)) continue;
351 
352  if (theMaxElement!=0 && result.size() >= theMaxElement) {
353  result.clear();
354  edm::LogError("TooManyTriplets")<<" number of triples exceed maximum. no triplets produced.";
355  return;
356  }
357  result.emplace_back( doublets.hit(ip,HitDoublets::inner), doublets.hit(ip,HitDoublets::outer), hits.theHits[KDdata].hit());
358  }
359  }
360  }
361  // std::cout << "found triplets " << result.size() << std::endl;
362 }
void swap(ora::Record &rh, ora::Record &lh)
Definition: Record.h:70
int i
Definition: DBlmapReader.cc:9
void build(std::vector< KDTreeNodeInfo > &eltList, const KDTreeBox &region)
T perp() const
Definition: PV3DBase.h:72
T y() const
Definition: PV3DBase.h:63
T r() const
Radius, same as mag()
bool ev
float fpi()
Definition: Pi.h:35
void search(const KDTreeBox &searchBox, std::vector< KDTreeNodeInfo > &resRecHitList)
U second(std::pair< T, U > const &p)
int seqNum() const
Definition: DetLayer.h:36
bool checkPhiInRange(float phi, float phi1, float phi2) const
T curvature(T InversePt, const edm::EventSetup &iSetup)
std::vector< SeedingLayerSetsHits::SeedingLayer > theLayers
T z() const
Definition: PV3DBase.h:64
tuple result
Definition: query.py:137
double f[11][100]
virtual HitDoublets doublets(const TrackingRegion &reg, const edm::Event &ev, const edm::EventSetup &es)
T min(T a, T b)
Definition: MathUtil.h:58
std::pair< float, float > mergePhiRanges(const std::pair< float, float > &r1, const std::pair< float, float > &r2) const
bool isBarrel() const
Definition: DetLayer.h:32
PixelRecoRange< float > Range
static double maxCurvature(const ThirdHitPredictionFromCircle *circle, double z1, double z2, double z3)
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:86
T x() const
Definition: PV3DBase.h:62
static const float fnSigmaRZ
DetLayer const * detLayer(layer l) const
tuple size
Write out results.
float ftwoPi()
Definition: Pi.h:36
double p3[4]
Definition: TauolaWrapper.h:91
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11
void PixelTripletLargeTipGenerator::init ( const HitPairGenerator pairs,
LayerCacheType layerCache 
)
overridevirtual

Implements HitTripletGeneratorFromPairAndLayers.

Definition at line 56 of file PixelTripletLargeTipGenerator.cc.

References HitPairGenerator::clone(), theLayerCache, and thePairGenerator.

58 {
59  thePairGenerator = pairs.clone();
60  theLayerCache = layerCache;
61 }
virtual HitPairGenerator * clone() const =0
std::pair< float, float > PixelTripletLargeTipGenerator::mergePhiRanges ( const std::pair< float, float > &  r1,
const std::pair< float, float > &  r2 
) const
private

Definition at line 371 of file PixelTripletLargeTipGenerator.cc.

References M_PI, bookConverter::max, and min().

Referenced by hitTriplets().

373 { float r2Min = r2.first;
374  float r2Max = r2.second;
375  while (r1.first - r2Min > +M_PI) r2Min += 2. * M_PI, r2Max += 2. * M_PI;
376  while (r1.first - r2Min < -M_PI) r2Min -= 2. * M_PI, r2Max -= 2. * M_PI;
377  return std::make_pair(min(r1.first, r2Min), max(r1.second, r2Max));
378 }
T min(T a, T b)
Definition: MathUtil.h:58
#define M_PI
const HitPairGenerator& PixelTripletLargeTipGenerator::pairGenerator ( ) const
inline

Definition at line 38 of file PixelTripletLargeTipGenerator.h.

References thePairGenerator.

38 { return *thePairGenerator; }
void PixelTripletLargeTipGenerator::setSeedingLayers ( SeedingLayerSetsHits::SeedingLayerSet  pairLayers,
std::vector< SeedingLayerSetsHits::SeedingLayer thirdLayers 
)
overridevirtual

Implements HitTripletGeneratorFromPairAndLayers.

Definition at line 63 of file PixelTripletLargeTipGenerator.cc.

References HitPairGenerator::setSeedingLayers(), theLayers, and thePairGenerator.

64  {
66  theLayers = thirdLayers;
67 }
std::vector< SeedingLayerSetsHits::SeedingLayer > theLayers
virtual void setSeedingLayers(SeedingLayerSetsHits::SeedingLayerSet layers)=0

Member Data Documentation

float PixelTripletLargeTipGenerator::dphi
private

Definition at line 57 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets(), and PixelTripletLargeTipGenerator().

float PixelTripletLargeTipGenerator::extraHitRPhitolerance
private

Definition at line 54 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

float PixelTripletLargeTipGenerator::extraHitRZtolerance
private

Definition at line 53 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

LayerCacheType* PixelTripletLargeTipGenerator::theLayerCache
private

Definition at line 50 of file PixelTripletLargeTipGenerator.h.

Referenced by init().

std::vector<SeedingLayerSetsHits::SeedingLayer> PixelTripletLargeTipGenerator::theLayers
private

Definition at line 49 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets(), and setSeedingLayers().

HitPairGenerator* PixelTripletLargeTipGenerator::thePairGenerator
private
bool PixelTripletLargeTipGenerator::useBend
private

Definition at line 56 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

bool PixelTripletLargeTipGenerator::useFixedPreFiltering
private

Definition at line 52 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets(), and PixelTripletLargeTipGenerator().

bool PixelTripletLargeTipGenerator::useMScat
private

Definition at line 55 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().