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 104 of file PixelTripletLargeTipGenerator.cc.

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

Referenced by hitTriplets().

109  {
110  auto const& doublets = thePairGenerator->doublets(region, ev, es, pairLayers);
111 
112  if (not doublets or doublets->empty())
113  return;
114 
116  hitTriplets(region, result, ev, es, *doublets, thirdLayers, nullptr, *theLayerCache);
117 }
assert(be >=bs)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
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 119 of file PixelTripletLargeTipGenerator.cc.

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

126  {
127  int size = thirdLayers.size();
128  const RecHitsSortedInPhi* thirdHitMap[size];
129  vector<const DetLayer*> thirdLayerDetLayer(size, nullptr);
130  for (int il = 0; il < size; ++il) {
131  thirdHitMap[il] = &layerCache(thirdLayers[il], region);
132  thirdLayerDetLayer[il] = thirdLayers[il].detLayer();
133  }
134  hitTriplets(region, result, es, doublets, thirdHitMap, thirdLayerDetLayer, size, tripletLastLayerIndex);
135 }
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 137 of file PixelTripletLargeTipGenerator.cc.

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

143  {
144  hitTriplets(region, result, es, doublets, thirdHitMap, thirdLayerDetLayer, nThirdLayers, nullptr);
145 }
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 147 of file PixelTripletLargeTipGenerator.cc.

References angle(), checkPhiInRange(), ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), pfMETCorrectionType0_cfi::correction, PixelRecoUtilities::curvature(), declareDynArray, HLT_2023v12_cff::doublets, dphi, relativeConstraints::empty, submitPVResolutionJobs::err, extraHitRPhitolerance, extraHitRZtolerance, f, fieldESToken_, dqmdumpme::first, nano_mu_digi_cff::float, fnSigmaRZ, Geom::ftwoPi(), edm::EventSetup::getData(), hfClusterShapes_cfi::hits, mps_fire::i, HitDoublets::inner, nano_mu_digi_cff::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, nano_mu_digi_cff::region, mps_fire::result, edm::second(), edm::swap(), HitTripletGeneratorFromPairAndLayers::theMaxElement, trackerTopologyESToken_, useBend, useFixedPreFiltering, useMScat, findQualityFiles::v, protons_cff::xi, and PV3DBase< T, PVType, FrameType >::z().

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