CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Types | Private Attributes
PixelTripletHLTGenerator Class Reference

#include <PixelTripletHLTGenerator.h>

Inheritance diagram for PixelTripletHLTGenerator:
HitTripletGeneratorFromPairAndLayers

Public Member Functions

void hitTriplets (const TrackingRegion &region, OrderedHitTriplets &trs, const edm::Event &ev, const edm::EventSetup &es, const SeedingLayerSetsHits::SeedingLayerSet &pairLayers, const std::vector< SeedingLayerSetsHits::SeedingLayer > &thirdLayers) override
 
void hitTriplets (const TrackingRegion &region, OrderedHitTriplets &trs, const edm::Event &ev, const edm::EventSetup &es, const HitDoublets &doublets, const std::vector< SeedingLayerSetsHits::SeedingLayer > &thirdLayers, std::vector< int > *tripletLastLayerIndex, LayerCacheType &layerCache)
 
void hitTriplets (const TrackingRegion &region, OrderedHitTriplets &result, const edm::EventSetup &es, const HitDoublets &doublets, const RecHitsSortedInPhi **thirdHitMap, const std::vector< const DetLayer * > &thirdLayerDetLayer, const int nThirdLayers) override
 
void hitTriplets (const TrackingRegion &region, OrderedHitTriplets &result, const edm::EventSetup &es, const HitDoublets &doublets, const RecHitsSortedInPhi **thirdHitMap, const std::vector< const DetLayer * > &thirdLayerDetLayer, const int nThirdLayers, std::vector< int > *tripletLastLayerIndex)
 
 PixelTripletHLTGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
 
 PixelTripletHLTGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
 
 ~PixelTripletHLTGenerator () override
 
- Public Member Functions inherited from HitTripletGeneratorFromPairAndLayers
 HitTripletGeneratorFromPairAndLayers (unsigned int maxElement=0)
 
 HitTripletGeneratorFromPairAndLayers (const edm::ParameterSet &pset)
 
void init (std::unique_ptr< HitPairGeneratorFromLayerPair > &&pairs, LayerCacheType *layerCache)
 
const HitPairGeneratorFromLayerPairpairGenerator () const
 
virtual ~HitTripletGeneratorFromPairAndLayers ()
 

Static Public Member Functions

static void fillDescriptions (edm::ParameterSetDescription &desc)
 
static const char * fillDescriptionsLabel ()
 
- Static Public Member Functions inherited from HitTripletGeneratorFromPairAndLayers
static void fillDescriptions (edm::ParameterSetDescription &desc)
 

Private Types

typedef CombinedHitTripletGenerator::LayerCacheType LayerCacheType
 

Private Attributes

const float dphi
 
const float extraHitRPhitolerance
 
const float extraHitRZtolerance
 
std::unique_ptr< SeedComparitortheComparitor
 
const bool useBend
 
const bool useFixedPreFiltering
 
const bool useMScat
 

Additional Inherited Members

- Public Types inherited from HitTripletGeneratorFromPairAndLayers
typedef LayerHitMapCache LayerCacheType
 
- Protected Attributes inherited from HitTripletGeneratorFromPairAndLayers
LayerCacheTypetheLayerCache
 
const unsigned int theMaxElement
 
std::unique_ptr< HitPairGeneratorFromLayerPairthePairGenerator
 

Detailed Description

Definition at line 20 of file PixelTripletHLTGenerator.h.

Member Typedef Documentation

Definition at line 22 of file PixelTripletHLTGenerator.h.

Constructor & Destructor Documentation

PixelTripletHLTGenerator::PixelTripletHLTGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector &&  iC 
)
inline

Definition at line 25 of file PixelTripletHLTGenerator.h.

References looper::cfg, fillDescriptions(), and ~PixelTripletHLTGenerator().

25 : PixelTripletHLTGenerator(cfg, iC) {}
PixelTripletHLTGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
PixelTripletHLTGenerator::PixelTripletHLTGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)

Definition at line 38 of file PixelTripletHLTGenerator.cc.

References beamerCreator::create(), reco::get(), edm::ParameterSet::getParameter(), AlCaHLTBitMon_QueryRunRegistry::string, and theComparitor.

40  useFixedPreFiltering(cfg.getParameter<bool>("useFixedPreFiltering")),
41  extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
42  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
43  useMScat(cfg.getParameter<bool>("useMultScattering")),
44  useBend(cfg.getParameter<bool>("useBending")),
45  dphi(useFixedPreFiltering ? cfg.getParameter<double>("phiPreFiltering") : 0)
46 {
47  edm::ParameterSet comparitorPSet =
48  cfg.getParameter<edm::ParameterSet>("SeedComparitorPSet");
49  std::string comparitorName = comparitorPSet.getParameter<std::string>("ComponentName");
50  if(comparitorName != "none") {
51  theComparitor.reset( SeedComparitorFactory::get()->create( comparitorName, comparitorPSet, iC) );
52  }
53 }
T getParameter(std::string const &) const
def create(alignables, pedeDump, additionalData, outputFile, config)
std::unique_ptr< SeedComparitor > theComparitor
T get(const Candidate &c)
Definition: component.h:55
PixelTripletHLTGenerator::~PixelTripletHLTGenerator ( )
override

Definition at line 55 of file PixelTripletHLTGenerator.cc.

Referenced by PixelTripletHLTGenerator().

55 {}

Member Function Documentation

void PixelTripletHLTGenerator::fillDescriptions ( edm::ParameterSetDescription desc)
static

Definition at line 57 of file PixelTripletHLTGenerator.cc.

References edm::ParameterSetDescription::add(), HitTripletGeneratorFromPairAndLayers::fillDescriptions(), edm::ParameterSetDescription::setAllowAnything(), and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by PixelTripletHLTGenerator().

57  {
59  desc.add<double>("extraHitRPhitolerance", 0.032);
60  desc.add<double>("extraHitRZtolerance", 0.037);
61  desc.add<bool>("useMultScattering", true);
62  desc.add<bool>("useBending", true);
63  desc.add<bool>("useFixedPreFiltering", false);
64  desc.add<double>("phiPreFiltering", 0.3);
65 
66  edm::ParameterSetDescription descComparitor;
67  descComparitor.add<std::string>("ComponentName", "none");
68  descComparitor.setAllowAnything(); // until we have moved SeedComparitor too to EDProducers
69  desc.add<edm::ParameterSetDescription>("SeedComparitorPSet", descComparitor);
70 }
void setAllowAnything()
allow any parameter label/value pairs
static void fillDescriptions(edm::ParameterSetDescription &desc)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static const char* PixelTripletHLTGenerator::fillDescriptionsLabel ( )
inlinestatic

Definition at line 31 of file PixelTripletHLTGenerator.h.

References HiRegitMuonDetachedTripletStep_cff::doublets, ev, hitTriplets(), and mps_fire::result.

31 { return "pixelTripletHLT"; }
void PixelTripletHLTGenerator::hitTriplets ( const TrackingRegion region,
OrderedHitTriplets trs,
const edm::Event ev,
const edm::EventSetup es,
const SeedingLayerSetsHits::SeedingLayerSet pairLayers,
const std::vector< SeedingLayerSetsHits::SeedingLayer > &  thirdLayers 
)
overridevirtual

Implements HitTripletGeneratorFromPairAndLayers.

Definition at line 72 of file PixelTripletHLTGenerator.cc.

References HiRegitMuonDetachedTripletStep_cff::doublets, HitTripletGeneratorFromPairAndLayers::theLayerCache, and HitTripletGeneratorFromPairAndLayers::thePairGenerator.

Referenced by fillDescriptionsLabel(), and hitTriplets().

78 {
79  auto const & doublets = thePairGenerator->doublets(region,ev,es, pairLayers);
80  if (doublets.empty()) return;
81 
82  assert(theLayerCache);
83  hitTriplets(region, result, ev, es, doublets, thirdLayers, nullptr, *theLayerCache);
84 }
void hitTriplets(const TrackingRegion &region, OrderedHitTriplets &trs, const edm::Event &ev, const edm::EventSetup &es, const SeedingLayerSetsHits::SeedingLayerSet &pairLayers, const std::vector< SeedingLayerSetsHits::SeedingLayer > &thirdLayers) override
std::unique_ptr< HitPairGeneratorFromLayerPair > thePairGenerator
void PixelTripletHLTGenerator::hitTriplets ( const TrackingRegion region,
OrderedHitTriplets trs,
const edm::Event ev,
const edm::EventSetup es,
const HitDoublets doublets,
const std::vector< SeedingLayerSetsHits::SeedingLayer > &  thirdLayers,
std::vector< int > *  tripletLastLayerIndex,
LayerCacheType layerCache 
)

Definition at line 86 of file PixelTripletHLTGenerator.cc.

References hitTriplets(), findQualityFiles::size, and theComparitor.

91  {
92  if (theComparitor) theComparitor->init(ev, es);
93 
94  int size = thirdLayers.size();
95  const RecHitsSortedInPhi * thirdHitMap[size];
96  vector<const DetLayer *> thirdLayerDetLayer(size,0);
97  for (int il=0; il<size; ++il)
98  {
99  thirdHitMap[il] = &layerCache(thirdLayers[il], region, es);
100  thirdLayerDetLayer[il] = thirdLayers[il].detLayer();
101  }
102  hitTriplets(region, result, es, doublets, thirdHitMap, thirdLayerDetLayer, size, tripletLastLayerIndex);
103 }
size
Write out results.
std::unique_ptr< SeedComparitor > theComparitor
void hitTriplets(const TrackingRegion &region, OrderedHitTriplets &trs, const edm::Event &ev, const edm::EventSetup &es, const SeedingLayerSetsHits::SeedingLayerSet &pairLayers, const std::vector< SeedingLayerSetsHits::SeedingLayer > &thirdLayers) override
void PixelTripletHLTGenerator::hitTriplets ( const TrackingRegion region,
OrderedHitTriplets result,
const edm::EventSetup es,
const HitDoublets doublets,
const RecHitsSortedInPhi **  thirdHitMap,
const std::vector< const DetLayer * > &  thirdLayerDetLayer,
const int  nThirdLayers 
)
overridevirtual

Implements HitTripletGeneratorFromPairAndLayers.

Definition at line 105 of file PixelTripletHLTGenerator.cc.

References hitTriplets().

113 {
114  hitTriplets(region, result, es, doublets, thirdHitMap, thirdLayerDetLayer, nThirdLayers, nullptr);
115 }
void hitTriplets(const TrackingRegion &region, OrderedHitTriplets &trs, const edm::Event &ev, const edm::EventSetup &es, const SeedingLayerSetsHits::SeedingLayerSet &pairLayers, const std::vector< SeedingLayerSetsHits::SeedingLayer > &thirdLayers) override
void PixelTripletHLTGenerator::hitTriplets ( const TrackingRegion region,
OrderedHitTriplets result,
const edm::EventSetup es,
const HitDoublets doublets,
const RecHitsSortedInPhi **  thirdHitMap,
const std::vector< const DetLayer * > &  thirdLayerDetLayer,
const int  nThirdLayers,
std::vector< int > *  tripletLastLayerIndex 
)

Definition at line 117 of file PixelTripletHLTGenerator.cc.

References angle(), checkPhiInRange(), constexpr, PixelRecoUtilities::curvature(), declareDynArray, HitDoublets::detLayer(), dphi, relativeConstraints::empty, extraHitRPhitolerance, extraHitRZtolerance, f, objects.autophobj::float, Geom::ftwoPi(), HitDoublets::hit(), hfClusterShapes_cfi::hits, mps_fire::i, ThirdHitRZPredictionBase::initLayer(), ThirdHitRZPredictionBase::initTolerance(), HitDoublets::inner, PixelRecoRange< T >::intersection(), DetLayer::isBarrel(), geometryCSVtoXML::line, LogDebug, M_PI, hpstanc_transforms::max, min(), normalizedPhi(), nSigmaPhi, nSigmaRZ, TrackingRegion::origin(), TrackingRegion::originRBound(), HitDoublets::outer, PV3DBase< T, PVType, FrameType >::perp(), HitDoublets::phi(), GeometricSearchDet::position(), proxim(), TrackingRegion::ptMin(), TCMET_cfi::radius, HitDoublets::rv(), DetLayer::seqNum(), OrderedHitTriplets::size(), HitDoublets::size(), edm::swap(), theComparitor, HitTripletGeneratorFromPairAndLayers::theMaxElement, useBend, useFixedPreFiltering, useMScat, findQualityFiles::v, PV3DBase< T, PVType, FrameType >::x(), HitDoublets::x(), hybridSuperClusters_cfi::xi, PV3DBase< T, PVType, FrameType >::y(), HitDoublets::y(), PV3DBase< T, PVType, FrameType >::z(), and HitDoublets::z().

123  {
124  auto outSeq = doublets.detLayer(HitDoublets::outer)->seqNum();
125 
126  float regOffset = region.origin().perp(); //try to take account of non-centrality (?)
127 
129  declareDynArray(ThirdHitCorrection, nThirdLayers, corrections);
130 
132 
133  using NodeInfo = KDTreeNodeInfo<unsigned int>;
134  std::vector<NodeInfo > layerTree; // re-used throughout
135  std::vector<unsigned int> foundNodes; // re-used thoughout
136  foundNodes.reserve(100);
137 
138  declareDynArray(KDTreeLinkerAlgo<unsigned int>,nThirdLayers, hitTree);
139  float rzError[nThirdLayers]; //save maximum errors
140 
141  const float maxDelphi = region.ptMin() < 0.3f ? float(M_PI)/4.f : float(M_PI)/8.f; // FIXME move to config??
142  const float maxphi = M_PI+maxDelphi, minphi = -maxphi; // increase to cater for any range
143  const float safePhi = M_PI-maxDelphi; // sideband
144 
145  // fill the prediction vector
146  for (int il=0; il<nThirdLayers; ++il) {
147  auto const & hits = *thirdHitMap[il];
148  ThirdHitRZPrediction<PixelRecoLineRZ> & pred = preds[il];
149  pred.initLayer(thirdLayerDetLayer[il]);
151 
152  corrections[il].init(es, region.ptMin(), *doublets.detLayer(HitDoublets::inner), *doublets.detLayer(HitDoublets::outer),
153  *thirdLayerDetLayer[il], useMScat, useBend);
154 
155  layerTree.clear();
156  float minv=999999.0f, maxv= -minv; // Initialise to extreme values in case no hits
157  float maxErr=0.0f;
158  for (unsigned int i=0; i!=hits.size(); ++i) {
159  auto angle = hits.phi(i);
160  auto v = hits.gv(i);
161  //use (phi,r) for endcaps rather than (phi,z)
162  minv = std::min(minv,v); maxv = std::max(maxv,v);
163  float myerr = hits.dv[i];
164  maxErr = std::max(maxErr,myerr);
165  layerTree.emplace_back(i, angle, v); // save it
166  // populate side-bands
167  if (angle>safePhi) layerTree.emplace_back(i, angle-Geom::ftwoPi(), v);
168  else if (angle<-safePhi) layerTree.emplace_back(i, angle+Geom::ftwoPi(), v);
169  }
170  KDTreeBox phiZ(minphi, maxphi, minv-0.01f, maxv+0.01f); // declare our bounds
171  //add fudge factors in case only one hit and also for floating-point inaccuracy
172  hitTree[il].build(layerTree, phiZ); // make KDtree
173  rzError[il] = maxErr; //save error
174  // std::cout << "layer " << thirdLayerDetLayer[il]->seqNum() << " " << layerTree.size() << std::endl;
175  }
176 
177  float imppar = region.originRBound();
178  float imppartmp = region.originRBound()+region.origin().perp();
179  float curv = PixelRecoUtilities::curvature(1.f/region.ptMin(), es);
180 
181  for (std::size_t ip =0; ip!=doublets.size(); ip++) {
182  auto xi = doublets.x(ip,HitDoublets::inner);
183  auto yi = doublets.y(ip,HitDoublets::inner);
184  auto zi = doublets.z(ip,HitDoublets::inner);
185  auto rvi = doublets.rv(ip,HitDoublets::inner);
186  auto xo = doublets.x(ip,HitDoublets::outer);
187  auto yo = doublets.y(ip,HitDoublets::outer);
188  auto zo = doublets.z(ip,HitDoublets::outer);
189  auto rvo = doublets.rv(ip,HitDoublets::outer);
190 
191  auto toPos = std::signbit(zo-zi);
192 
193  PixelRecoPointRZ point1(rvi, zi);
194  PixelRecoPointRZ point2(rvo, zo);
195  PixelRecoLineRZ line(point1, point2);
196  ThirdHitPredictionFromInvParabola predictionRPhi(xi-region.origin().x(),yi-region.origin().y(),
197  xo-region.origin().x(),yo-region.origin().y(),
198  imppar,curv,extraHitRPhitolerance);
199 
200  ThirdHitPredictionFromInvParabola predictionRPhitmp(xi,yi,xo,yo,imppartmp,curv,extraHitRPhitolerance);
201 
202  // printf("++Constr %f %f %f %f %f %f %f\n",xi,yi,xo,yo,imppartmp,curv,extraHitRPhitolerance);
203 
204  // std::cout << ip << ": " << point1.r() << ","<< point1.z() << " "
205  // << point2.r() << ","<< point2.z() <<std::endl;
206 
207  for (int il=0; il!=nThirdLayers; ++il) {
208  const DetLayer * layer = thirdLayerDetLayer[il];
209  auto barrelLayer = layer->isBarrel();
210 
211  if ( (!barrelLayer) & (toPos != std::signbit(layer->position().z())) ) continue;
212 
213  if (hitTree[il].empty()) continue; // Don't bother if no hits
214 
215  auto const & hits = *thirdHitMap[il];
216 
217  auto & correction = corrections[il];
218 
219  correction.init(line, point2, outSeq);
220 
221  auto & predictionRZ = preds[il];
222 
223  predictionRZ.initPropagator(&line);
224  Range rzRange = predictionRZ();
225  correction.correctRZRange(rzRange);
226 
227  Range phiRange;
228  if (useFixedPreFiltering) {
229  float phi0 = doublets.phi(ip,HitDoublets::outer);
230  phiRange = Range(phi0-dphi,phi0+dphi);
231  }
232  else {
233  Range radius;
234  if (barrelLayer) {
235  radius = predictionRZ.detRange();
236  } else {
237  radius = Range(max(rzRange.min(), predictionRZ.detSize().min()),
238  min(rzRange.max(), predictionRZ.detSize().max()) );
239  }
240  if (radius.empty()) continue;
241 
242  // std::cout << "++R " << radius.min() << " " << radius.max() << std::endl;
243 
244  auto rPhi1 = predictionRPhitmp(radius.max());
245  bool ok1 = !rPhi1.empty();
246  if (ok1) {
247  correction.correctRPhiRange(rPhi1);
248  rPhi1.first /= radius.max();
249  rPhi1.second /= radius.max();
250  }
251  auto rPhi2 = predictionRPhitmp(radius.min());
252  bool ok2 = !rPhi2.empty();
253  if (ok2) {
254  correction.correctRPhiRange(rPhi2);
255  rPhi2.first /= radius.min();
256  rPhi2.second /= radius.min();
257  }
258 
259  if (ok1) {
260  rPhi1.first = normalizedPhi(rPhi1.first);
261  rPhi1.second = proxim(rPhi1.second,rPhi1.first);
262  if(ok2) {
263  rPhi2.first = proxim(rPhi2.first,rPhi1.first);
264  rPhi2.second = proxim(rPhi2.second,rPhi1.first);
265  phiRange = rPhi1.sum(rPhi2);
266  } else phiRange=rPhi1;
267  } else if(ok2) {
268  rPhi2.first = normalizedPhi(rPhi2.first);
269  rPhi2.second = proxim(rPhi2.second,rPhi2.first);
270  phiRange=rPhi2;
271  } else continue;
272 
273  }
274 
275  constexpr float nSigmaRZ = 3.46410161514f; // std::sqrt(12.f); // ...and continue as before
276  constexpr float nSigmaPhi = 3.f;
277 
278  foundNodes.clear(); // Now recover hits in bounding box...
279  float prmin=phiRange.min(), prmax=phiRange.max();
280 
281 
282  if (prmax-prmin>maxDelphi) {
283  auto prm = phiRange.mean();
284  prmin = prm - 0.5f*maxDelphi;
285  prmax = prm + 0.5f*maxDelphi;
286  }
287 
288  if (barrelLayer)
289  {
290  Range regMax = predictionRZ.detRange();
291  Range regMin = predictionRZ(regMax.min()-regOffset);
292  regMax = predictionRZ(regMax.max()+regOffset);
293  correction.correctRZRange(regMin);
294  correction.correctRZRange(regMax);
295  if (regMax.min() < regMin.min()) { swap(regMax, regMin);}
296  KDTreeBox phiZ(prmin, prmax, regMin.min()-nSigmaRZ*rzError[il], regMax.max()+nSigmaRZ*rzError[il]);
297  hitTree[il].search(phiZ, foundNodes);
298  }
299  else
300  {
301  KDTreeBox phiZ(prmin, prmax,
302  rzRange.min()-regOffset-nSigmaRZ*rzError[il],
303  rzRange.max()+regOffset+nSigmaRZ*rzError[il]);
304  hitTree[il].search(phiZ, foundNodes);
305  }
306 
307  // std::cout << ip << ": " << thirdLayerDetLayer[il]->seqNum() << " " << foundNodes.size() << " " << prmin << " " << prmax << std::endl;
308 
309 
310  // int kk=0;
311  for (auto KDdata : foundNodes) {
312 
313  if (theMaxElement!=0 && result.size() >= theMaxElement){
314  result.clear();
315  if(tripletLastLayerIndex) tripletLastLayerIndex->clear();
316  edm::LogError("TooManyTriplets")<<" number of triples exceeds maximum. no triplets produced.";
317  return;
318  }
319 
320  float p3_u = hits.u[KDdata];
321  float p3_v = hits.v[KDdata];
322  float p3_phi = hits.lphi[KDdata];
323 
324  //if ((kk++)%100==0)
325  //std::cout << kk << ": " << p3_u << " " << p3_v << " " << p3_phi << std::endl;
326 
327 
328  Range allowed = predictionRZ(p3_u);
329  correction.correctRZRange(allowed);
330  float vErr = nSigmaRZ *hits.dv[KDdata];
331  Range hitRange(p3_v-vErr, p3_v+vErr);
332  Range crossingRange = allowed.intersection(hitRange);
333  if (crossingRange.empty()) continue;
334 
335  float ir = 1.f/hits.rv(KDdata);
336  // limit error to 90 degree
337  constexpr float maxPhiErr = 0.5*M_PI;
338  float phiErr = nSigmaPhi * hits.drphi[KDdata]*ir;
339  phiErr = std::min(maxPhiErr, phiErr);
340  bool nook=true;
341  for (int icharge=-1; icharge <=1; icharge+=2) {
342  Range rangeRPhi = predictionRPhi(hits.rv(KDdata), icharge);
343  if(rangeRPhi.first>rangeRPhi.second) continue; // range is empty
344  correction.correctRPhiRange(rangeRPhi);
345  if (checkPhiInRange(p3_phi, rangeRPhi.first*ir-phiErr, rangeRPhi.second*ir+phiErr,maxPhiErr)) {
346  // insert here check with comparitor
347  OrderedHitTriplet hittriplet( doublets.hit(ip,HitDoublets::inner), doublets.hit(ip,HitDoublets::outer), hits.theHits[KDdata].hit());
348  if (!theComparitor || theComparitor->compatible(hittriplet) ) {
349  result.push_back( hittriplet );
350  // to bookkeep the triplets and 3rd layers in triplet EDProducer
351  if(tripletLastLayerIndex) tripletLastLayerIndex->push_back(il);
352  } else {
353  LogDebug("RejectedTriplet") << "rejected triplet from comparitor ";
354  }
355  nook=false; break;
356  }
357  }
358  if (nook) LogDebug("RejectedTriplet") << "rejected triplet from second phicheck " << p3_phi;
359  }
360  }
361  }
362  // std::cout << "triplets " << result.size() << std::endl;
363 }
#define LogDebug(id)
float originRBound() const
bounds the particle vertex in the transverse plane
std::size_t size() const
T perp() const
Definition: PV3DBase.h:72
GlobalPoint const & origin() const
constexpr float ftwoPi()
Definition: Pi.h:36
void initLayer(const DetLayer *layer)
PixelRecoRange< float > Range
float x(int i, layer l) const
T y() const
Definition: PV3DBase.h:63
T normalizedPhi(T phi)
Definition: normalizedPhi.h:9
#define constexpr
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:116
std::unique_ptr< SeedComparitor > theComparitor
float z(int i, layer l) const
int seqNum() const
Definition: DetLayer.h:36
T curvature(T InversePt, const edm::EventSetup &iSetup)
bool checkPhiInRange(T phi, T phi1, T phi2, float maxDphi=float(M_PI))
Definition: normalizedPhi.h:34
T z() const
Definition: PV3DBase.h:64
double f[11][100]
T min(T a, T b)
Definition: MathUtil.h:58
float y(int i, layer l) const
unsigned int size() const override
SeedingHitSet::ConstRecHitPointer Hit
BaseTrackerRecHit const * Hit
#define M_PI
bool isBarrel() const
Definition: DetLayer.h:32
float phi(int i, layer l) const
virtual const Surface::PositionType & position() const
Returns position of the surface.
float ptMin() const
minimal pt of interest
T proxim(T b, T a)
Definition: normalizedPhi.h:14
Hit const & hit(int i, layer l) const
PixelRecoRange< T > intersection(const PixelRecoRange< T > &r) const
#define declareDynArray(T, n, x)
Definition: DynArray.h:59
void initTolerance(float tolerance)
T x() const
Definition: PV3DBase.h:62
DetLayer const * detLayer(layer l) const
float rv(int i, layer l) const
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11

Member Data Documentation

const float PixelTripletHLTGenerator::dphi
private

Definition at line 68 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

const float PixelTripletHLTGenerator::extraHitRPhitolerance
private

Definition at line 65 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

const float PixelTripletHLTGenerator::extraHitRZtolerance
private

Definition at line 64 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

std::unique_ptr<SeedComparitor> PixelTripletHLTGenerator::theComparitor
private

Definition at line 69 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets(), and PixelTripletHLTGenerator().

const bool PixelTripletHLTGenerator::useBend
private

Definition at line 67 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

const bool PixelTripletHLTGenerator::useFixedPreFiltering
private

Definition at line 63 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

const bool PixelTripletHLTGenerator::useMScat
private

Definition at line 66 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().