CMS 3D CMS Logo

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

#include <PixelTripletLargeTipGenerator.h>

Inheritance diagram for PixelTripletLargeTipGenerator:
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)
 
 PixelTripletLargeTipGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
 
 PixelTripletLargeTipGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
 
 ~PixelTripletLargeTipGenerator () 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
 
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecordfieldESToken_
 
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecordmsmakerESToken_
 
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcdtrackerTopologyESToken_
 
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

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 23 of file PixelTripletLargeTipGenerator.h.

Member Typedef Documentation

◆ LayerCacheType

Definition at line 24 of file PixelTripletLargeTipGenerator.h.

Constructor & Destructor Documentation

◆ PixelTripletLargeTipGenerator() [1/2]

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

Definition at line 27 of file PixelTripletLargeTipGenerator.h.

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

◆ PixelTripletLargeTipGenerator() [2/2]

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

Definition at line 47 of file PixelTripletLargeTipGenerator.cc.

References edm::ConsumesCollector::esConsumes(), msmakerESToken_, and useMScat.

49  useFixedPreFiltering(cfg.getParameter<bool>("useFixedPreFiltering")),
50  extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
51  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
52  useMScat(cfg.getParameter<bool>("useMultScattering")),
53  useBend(cfg.getParameter<bool>("useBending")),
54  dphi(useFixedPreFiltering ? cfg.getParameter<double>("phiPreFiltering") : 0),
57  if (useMScat) {
59  }
60 }
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > fieldESToken_
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > trackerTopologyESToken_
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > msmakerESToken_

◆ ~PixelTripletLargeTipGenerator()

PixelTripletLargeTipGenerator::~PixelTripletLargeTipGenerator ( )
override

Definition at line 62 of file PixelTripletLargeTipGenerator.cc.

62 {}

Member Function Documentation

◆ fillDescriptions()

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

Definition at line 64 of file PixelTripletLargeTipGenerator.cc.

References submitPVResolutionJobs::desc, and HitTripletGeneratorFromPairAndLayers::fillDescriptions().

64  {
66  // Defaults for the extraHit*tolerance are set to 0 here since that
67  // was already the case in practice in all offline occurrances.
68  desc.add<double>("extraHitRPhitolerance", 0); // default in old python was 0.032
69  desc.add<double>("extraHitRZtolerance", 0); // default in old python was 0.037
70  desc.add<bool>("useMultScattering", true);
71  desc.add<bool>("useBending", true);
72  desc.add<bool>("useFixedPreFiltering", false);
73  desc.add<double>("phiPreFiltering", 0.3);
74 }
static void fillDescriptions(edm::ParameterSetDescription &desc)

◆ fillDescriptionsLabel()

static const char* PixelTripletLargeTipGenerator::fillDescriptionsLabel ( )
inlinestatic

Definition at line 34 of file PixelTripletLargeTipGenerator.h.

34 { return "pixelTripletLargeTip"; }

◆ hitTriplets() [1/4]

void PixelTripletLargeTipGenerator::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 88 of file PixelTripletLargeTipGenerator.cc.

References cms::cuda::assert(), HLT_2023v12_cff::doublets, makeMEIFBenchmarkPlots::ev, HLT_2023v12_cff::region, mps_fire::result, HitTripletGeneratorFromPairAndLayers::theLayerCache, and HitTripletGeneratorFromPairAndLayers::thePairGenerator.

Referenced by hitTriplets().

93  {
94  auto const& doublets = thePairGenerator->doublets(region, ev, es, pairLayers);
95 
96  if (doublets.empty())
97  return;
98 
100  hitTriplets(region, result, ev, es, doublets, thirdLayers, nullptr, *theLayerCache);
101 }
assert(be >=bs)
std::unique_ptr< HitPairGeneratorFromLayerPair > thePairGenerator
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

◆ hitTriplets() [2/4]

void PixelTripletLargeTipGenerator::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 103 of file PixelTripletLargeTipGenerator.cc.

References HLT_2023v12_cff::doublets, hitTriplets(), HLT_2023v12_cff::region, mps_fire::result, and findQualityFiles::size.

110  {
111  int size = thirdLayers.size();
112  const RecHitsSortedInPhi* thirdHitMap[size];
113  vector<const DetLayer*> thirdLayerDetLayer(size, nullptr);
114  for (int il = 0; il < size; ++il) {
115  thirdHitMap[il] = &layerCache(thirdLayers[il], region);
116  thirdLayerDetLayer[il] = thirdLayers[il].detLayer();
117  }
118  hitTriplets(region, result, es, doublets, thirdHitMap, thirdLayerDetLayer, size, tripletLastLayerIndex);
119 }
size
Write out results.
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

◆ hitTriplets() [3/4]

void PixelTripletLargeTipGenerator::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 121 of file PixelTripletLargeTipGenerator.cc.

References HLT_2023v12_cff::doublets, hitTriplets(), HLT_2023v12_cff::region, and mps_fire::result.

127  {
128  hitTriplets(region, result, es, doublets, thirdHitMap, thirdLayerDetLayer, nThirdLayers, nullptr);
129 }
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

◆ hitTriplets() [4/4]

void PixelTripletLargeTipGenerator::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 131 of file PixelTripletLargeTipGenerator.cc.

References angle(), checkPhiInRange(), pfMETCorrectionType0_cfi::correction, PixelRecoUtilities::curvature(), declareDynArray, HLT_2023v12_cff::doublets, dphi, relativeConstraints::empty, submitPVResolutionJobs::err, extraHitRPhitolerance, extraHitRZtolerance, f, fieldESToken_, dqmdumpme::first, dqmMemoryStats::float, fnSigmaRZ, Geom::ftwoPi(), edm::EventSetup::getData(), hfClusterShapes_cfi::hits, mps_fire::i, HitDoublets::inner, pixelTopology::layer, mps_splice::line, M_PI, SiStripPI::max, ThirdHitPredictionFromCircle::HelixRZ::maxCurvature(), SiStripPI::min, msmakerESToken_, normalizedPhi(), nSigmaPhi, nSigmaRZ, HitDoublets::outer, chargedHadronTrackResolutionFilter_cfi::p3, PV3DBase< T, PVType, FrameType >::perp(), proxim(), Basic2DVector< T >::r(), CosmicsPD_Skims::radius, HLT_2023v12_cff::region, mps_fire::result, edm::second(), std::swap(), HitTripletGeneratorFromPairAndLayers::theMaxElement, trackerTopologyESToken_, useBend, useFixedPreFiltering, useMScat, findQualityFiles::v, protons_cff::xi, and PV3DBase< T, PVType, FrameType >::z().

138  {
139  const TrackerTopology* tTopo = &es.getData(trackerTopologyESToken_);
140  const auto& field = es.getData(fieldESToken_);
141  const MultipleScatteringParametrisationMaker* msmaker = nullptr;
142  if (useMScat) {
143  msmaker = &es.getData(msmakerESToken_);
144  }
145 
146  auto outSeq = doublets.detLayer(HitDoublets::outer)->seqNum();
147 
148  using NodeInfo = KDTreeNodeInfo<unsigned int, 2>;
149  std::vector<NodeInfo> layerTree; // re-used throughout
150  std::vector<unsigned int> foundNodes; // re-used throughout
151  foundNodes.reserve(100);
152 
153  declareDynArray(KDTreeLinkerAlgo<unsigned int>, nThirdLayers, hitTree);
154  declareDynArray(LayerRZPredictions, nThirdLayers, mapPred);
155 
156  float rzError[nThirdLayers]; //save maximum errors
157 
158  const float maxDelphi = region.ptMin() < 0.3f ? float(M_PI) / 4.f : float(M_PI) / 8.f; // FIXME move to config??
159  const float maxphi = M_PI + maxDelphi, minphi = -maxphi; // increase to cater for any range
160  const float safePhi = M_PI - maxDelphi; // sideband
161 
162  for (int il = 0; il < nThirdLayers; il++) {
163  auto const& hits = *thirdHitMap[il];
164 
165  const DetLayer* layer = thirdLayerDetLayer[il];
166  LayerRZPredictions& predRZ = mapPred[il];
167  predRZ.line.initLayer(layer);
168  predRZ.helix1.initLayer(layer);
169  predRZ.helix2.initLayer(layer);
170  predRZ.line.initTolerance(extraHitRZtolerance);
171  predRZ.helix1.initTolerance(extraHitRZtolerance);
172  predRZ.helix2.initTolerance(extraHitRZtolerance);
173  predRZ.rzPositionFixup = MatchedHitRZCorrectionFromBending(layer, tTopo);
174  predRZ.correction.init(region.ptMin(),
175  *doublets.detLayer(HitDoublets::inner),
176  *doublets.detLayer(HitDoublets::outer),
177  *thirdLayerDetLayer[il],
178  useMScat,
179  msmaker,
180  false,
181  nullptr);
182 
183  layerTree.clear();
184  float minv = 999999.0;
185  float maxv = -999999.0; // Initialise to extreme values in case no hits
186  float maxErr = 0.0f;
187  for (unsigned int i = 0; i != hits.size(); ++i) {
188  auto angle = hits.phi(i);
189  auto v = hits.gv(i);
190  //use (phi,r) for endcaps rather than (phi,z)
191  minv = std::min(minv, v);
192  maxv = std::max(maxv, v);
193  float myerr = hits.dv[i];
194  maxErr = std::max(maxErr, myerr);
195  layerTree.emplace_back(i, angle, v); // save it
196  // populate side-bands
197  if (angle > safePhi)
198  layerTree.emplace_back(i, angle - Geom::ftwoPi(), v);
199  else if (angle < -safePhi)
200  layerTree.emplace_back(i, angle + Geom::ftwoPi(), v);
201  }
202  KDTreeBox phiZ(minphi, maxphi, minv - 0.01f, maxv + 0.01f); // declare our bounds
203  //add fudge factors in case only one hit and also for floating-point inaccuracy
204  hitTree[il].build(layerTree, phiZ); // make KDtree
205  rzError[il] = maxErr; //save error
206  }
207 
208  double curv = PixelRecoUtilities::curvature(1. / region.ptMin(), field);
209 
210  for (std::size_t ip = 0; ip != doublets.size(); ip++) {
211  auto xi = doublets.x(ip, HitDoublets::inner);
212  auto yi = doublets.y(ip, HitDoublets::inner);
213  auto zi = doublets.z(ip, HitDoublets::inner);
214  // auto rvi = doublets.rv(ip,HitDoublets::inner);
215  auto xo = doublets.x(ip, HitDoublets::outer);
216  auto yo = doublets.y(ip, HitDoublets::outer);
217  auto zo = doublets.z(ip, HitDoublets::outer);
218  // auto rvo = doublets.rv(ip,HitDoublets::outer);
219  GlobalPoint gp1(xi, yi, zi);
220  GlobalPoint gp2(xo, yo, zo);
221 
222  auto toPos = std::signbit(zo - zi);
223 
224  PixelRecoLineRZ line(gp1, gp2);
225  PixelRecoPointRZ point2(gp2.perp(), zo);
226  ThirdHitPredictionFromCircle predictionRPhi(gp1, gp2, extraHitRPhitolerance);
227 
228  Range generalCurvature = predictionRPhi.curvature(region.originRBound());
229  if (!intersect(generalCurvature, Range(-curv, curv)))
230  continue;
231 
232  for (int il = 0; il < nThirdLayers; il++) {
233  if (hitTree[il].empty())
234  continue; // Don't bother if no hits
235  const DetLayer* layer = thirdLayerDetLayer[il];
236  bool barrelLayer = layer->isBarrel();
237 
238  if ((!barrelLayer) & (toPos != std::signbit(layer->position().z())))
239  continue;
240 
241  Range curvature = generalCurvature;
242 
243  LayerRZPredictions& predRZ = mapPred[il];
244  predRZ.line.initPropagator(&line);
245 
246  auto& correction = predRZ.correction;
247  correction.init(line, point2, outSeq);
248 
249  Range rzRange;
250  if (useBend) {
251  // For the barrel region:
252  // swiping the helix passing through the two points across from
253  // negative to positive bending, can give us a sort of U-shaped
254  // projection onto the phi-z (barrel) or r-z plane (forward)
255  // so we checking minimum/maximum of all three possible extrema
256  //
257  // For the endcap region:
258  // Checking minimum/maximum radius of the helix projection
259  // onto an endcap plane, here we have to guard against
260  // looping tracks, when phi(delta z) gets out of control.
261  // HelixRZ::rAtZ should not follow looping tracks, but clamp
262  // to the minimum reachable r with the next-best lower |curvature|.
263  // So same procedure as for the barrel region can be applied.
264  //
265  // In order to avoid looking for potential looping tracks at all
266  // we also clamp the allowed curvature range for this layer,
267  // and potentially fail the layer entirely
268 
269  if (!barrelLayer) {
270  Range z3s = predRZ.line.detRange();
271  double z3 = z3s.first < 0 ? std::max(z3s.first, z3s.second) : std::min(z3s.first, z3s.second);
272  double maxCurvature = HelixRZ::maxCurvature(&predictionRPhi, gp1.z(), gp2.z(), z3);
273  if (!intersect(curvature, Range(-maxCurvature, maxCurvature)))
274  continue;
275  }
276 
277  HelixRZ helix1(&predictionRPhi, gp1.z(), gp2.z(), curvature.first);
278  HelixRZ helix2(&predictionRPhi, gp1.z(), gp2.z(), curvature.second);
279 
280  predRZ.helix1.initPropagator(&helix1);
281  predRZ.helix2.initPropagator(&helix2);
282 
283  Range rzRanges[2] = {predRZ.helix1(), predRZ.helix2()};
284  predRZ.helix1.initPropagator(nullptr);
285  predRZ.helix2.initPropagator(nullptr);
286 
287  rzRange.first = std::min(rzRanges[0].first, rzRanges[1].first);
288  rzRange.second = std::max(rzRanges[0].second, rzRanges[1].second);
289 
290  // if the allowed curvatures include a straight line,
291  // this can give us another extremum for allowed r/z
292  if (curvature.first * curvature.second < 0.0) {
293  Range rzLineRange = predRZ.line();
294  rzRange.first = std::min(rzRange.first, rzLineRange.first);
295  rzRange.second = std::max(rzRange.second, rzLineRange.second);
296  }
297  } else {
298  rzRange = predRZ.line();
299  }
300 
301  if (rzRange.first >= rzRange.second)
302  continue;
303 
304  correction.correctRZRange(rzRange);
305 
306  Range phiRange;
307  if (useFixedPreFiltering) {
308  float phi0 = doublets.phi(ip, HitDoublets::outer);
309  phiRange = Range(phi0 - dphi, phi0 + dphi);
310  } else {
311  Range radius;
312 
313  if (barrelLayer) {
314  radius = predRZ.line.detRange();
315  if (!intersect(rzRange, predRZ.line.detSize()))
316  continue;
317  } else {
318  radius = rzRange;
319  if (!intersect(radius, predRZ.line.detSize()))
320  continue;
321  }
322 
323  //gc: predictionRPhi uses the cosine rule to find the phi of the 3rd point at radius, assuming the pairCurvature range [-c,+c]
324  if ((curvature.first < 0.0f) & (curvature.second < 0.0f)) {
325  radius.swap();
326  } else if ((curvature.first >= 0.0f) & (curvature.second >= 0.0f)) {
327  ;
328  } else {
329  radius.first = radius.second;
330  }
331  auto phi12 = predictionRPhi.phi(curvature.first, radius.first);
332  auto phi22 = predictionRPhi.phi(curvature.second, radius.second);
333  phi12 = normalizedPhi(phi12);
334  phi22 = proxim(phi22, phi12);
335  phiRange = Range(phi12, phi22);
336  phiRange.sort();
337  auto rmean = radius.mean();
338  phiRange.first *= rmean;
339  phiRange.second *= rmean;
340  correction.correctRPhiRange(phiRange);
341  phiRange.first /= rmean;
342  phiRange.second /= rmean;
343  }
344 
345  foundNodes.clear(); // Now recover hits in bounding box...
346  float prmin = phiRange.min(), prmax = phiRange.max(); //get contiguous range
347 
348  if (prmax - prmin > maxDelphi) {
349  auto prm = phiRange.mean();
350  prmin = prm - 0.5f * maxDelphi;
351  prmax = prm + 0.5f * maxDelphi;
352  }
353 
354  if (barrelLayer) {
355  Range regMax = predRZ.line.detRange();
356  Range regMin = predRZ.line(regMax.min());
357  regMax = predRZ.line(regMax.max());
358  correction.correctRZRange(regMin);
359  correction.correctRZRange(regMax);
360  if (regMax.min() < regMin.min()) {
361  std::swap(regMax, regMin);
362  }
363  KDTreeBox phiZ(prmin, prmax, regMin.min() - fnSigmaRZ * rzError[il], regMax.max() + fnSigmaRZ * rzError[il]);
364  hitTree[il].search(phiZ, foundNodes);
365  } else {
366  KDTreeBox phiZ(prmin, prmax, rzRange.min() - fnSigmaRZ * rzError[il], rzRange.max() + fnSigmaRZ * rzError[il]);
367  hitTree[il].search(phiZ, foundNodes);
368  }
369 
370  MatchedHitRZCorrectionFromBending l2rzFixup(doublets.hit(ip, HitDoublets::outer)->det()->geographicalId(), tTopo);
371  MatchedHitRZCorrectionFromBending l3rzFixup = predRZ.rzPositionFixup;
372 
373  auto const& hits = *thirdHitMap[il];
374  for (auto KDdata : foundNodes) {
375  GlobalPoint p3 = hits.gp(KDdata);
376  double p3_r = p3.perp();
377  double p3_z = p3.z();
378  float p3_phi = hits.phi(KDdata);
379 
380  Range rangeRPhi = predictionRPhi(curvature, p3_r);
381  correction.correctRPhiRange(rangeRPhi);
382 
383  float ir = 1.f / p3_r;
384  // limit error to 90 degree
385  constexpr float maxPhiErr = 0.5 * M_PI;
386  float phiErr = nSigmaPhi * hits.drphi[KDdata] * ir;
387  phiErr = std::min(maxPhiErr, phiErr);
388  if (!checkPhiInRange(p3_phi, rangeRPhi.first * ir - phiErr, rangeRPhi.second * ir + phiErr, maxPhiErr))
389  continue;
390 
391  Basic2DVector<double> thc(p3.x(), p3.y());
392 
393  auto curv_ = predictionRPhi.curvature(thc);
394  double p2_r = point2.r();
395  double p2_z = point2.z(); // they will be modified!
396 
397  l2rzFixup(predictionRPhi, curv_, *doublets.hit(ip, HitDoublets::outer), p2_r, p2_z, tTopo);
398  l3rzFixup(predictionRPhi, curv_, *hits.theHits[KDdata].hit(), p3_r, p3_z, tTopo);
399 
400  Range rangeRZ;
401  if (useBend) {
402  HelixRZ updatedHelix(&predictionRPhi, gp1.z(), p2_z, curv_);
403  rangeRZ = predRZ.helix1(barrelLayer ? p3_r : p3_z, updatedHelix);
404  } else {
405  float tIP = predictionRPhi.transverseIP(thc);
406  PixelRecoPointRZ updatedPoint2(p2_r, p2_z);
407  PixelRecoLineRZ updatedLine(line.origin(), point2, tIP);
408  rangeRZ = predRZ.line(barrelLayer ? p3_r : p3_z, line);
409  }
410  correction.correctRZRange(rangeRZ);
411 
412  double err = nSigmaRZ * hits.dv[KDdata];
413 
414  rangeRZ.first -= err, rangeRZ.second += err;
415 
416  if (!rangeRZ.inside(barrelLayer ? p3_z : p3_r))
417  continue;
418 
419  if (theMaxElement != 0 && result.size() >= theMaxElement) {
420  result.clear();
421  if (tripletLastLayerIndex)
422  tripletLastLayerIndex->clear();
423  edm::LogError("TooManyTriplets") << " number of triples exceed maximum. no triplets produced.";
424  return;
425  }
426  result.emplace_back(
427  doublets.hit(ip, HitDoublets::inner), doublets.hit(ip, HitDoublets::outer), hits.theHits[KDdata].hit());
428  // to bookkeep the triplets and 3rd layers in triplet EDProducer
429  if (tripletLastLayerIndex)
430  tripletLastLayerIndex->push_back(il);
431  }
432  }
433  }
434  // std::cout << "found triplets " << result.size() << std::endl;
435 }
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
constexpr float ftwoPi()
Definition: Pi.h:36
constexpr T normalizedPhi(T phi)
Definition: normalizedPhi.h:8
Log< level::Error, false > LogError
T curvature(T InversePt, const MagneticField &field)
constexpr T proxim(T b, T a)
Definition: normalizedPhi.h:14
constexpr std::array< uint8_t, layerIndexSize< TrackerTraits > > layer
U second(std::pair< T, U > const &p)
constexpr float fnSigmaRZ
constexpr double nSigmaRZ
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
double f[11][100]
constexpr float nSigmaPhi
#define M_PI
T r() const
Radius, same as mag()
static double maxCurvature(const ThirdHitPredictionFromCircle *circle, double z1, double z2, double z3)
PixelRecoRange< float > Range
constexpr bool checkPhiInRange(T phi, T phi1, T phi2, float maxDphi=float(M_PI))
Definition: normalizedPhi.h:33
#define declareDynArray(T, n, x)
Definition: DynArray.h:91
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > fieldESToken_
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > trackerTopologyESToken_
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > msmakerESToken_
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11

Member Data Documentation

◆ dphi

const float PixelTripletLargeTipGenerator::dphi
private

Definition at line 75 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

◆ extraHitRPhitolerance

const float PixelTripletLargeTipGenerator::extraHitRPhitolerance
private

Definition at line 72 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

◆ extraHitRZtolerance

const float PixelTripletLargeTipGenerator::extraHitRZtolerance
private

Definition at line 71 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

◆ fieldESToken_

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

Definition at line 78 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

◆ msmakerESToken_

edm::ESGetToken<MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord> PixelTripletLargeTipGenerator::msmakerESToken_
private

Definition at line 79 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets(), and PixelTripletLargeTipGenerator().

◆ trackerTopologyESToken_

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

Definition at line 77 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

◆ useBend

const bool PixelTripletLargeTipGenerator::useBend
private

Definition at line 74 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

◆ useFixedPreFiltering

const bool PixelTripletLargeTipGenerator::useFixedPreFiltering
private

Definition at line 70 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

◆ useMScat

const bool PixelTripletLargeTipGenerator::useMScat
private

Definition at line 73 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets(), and PixelTripletLargeTipGenerator().