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
 
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecordfieldToken_
 
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecordmsmakerToken_
 
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 23 of file PixelTripletHLTGenerator.h.

Member Typedef Documentation

◆ LayerCacheType

Definition at line 24 of file PixelTripletHLTGenerator.h.

Constructor & Destructor Documentation

◆ PixelTripletHLTGenerator() [1/2]

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

Definition at line 27 of file PixelTripletHLTGenerator.h.

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

◆ PixelTripletHLTGenerator() [2/2]

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

Definition at line 32 of file PixelTripletHLTGenerator.cc.

References looper::cfg, edm::ConsumesCollector::esConsumes(), get, edm::ParameterSet::getParameter(), msmakerToken_, AlCaHLTBitMon_QueryRunRegistry::string, theComparitor, and useMScat.

34  fieldToken_(iC.esConsumes()),
35  useFixedPreFiltering(cfg.getParameter<bool>("useFixedPreFiltering")),
36  extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
37  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
38  useMScat(cfg.getParameter<bool>("useMultScattering")),
39  useBend(cfg.getParameter<bool>("useBending")),
40  dphi(useFixedPreFiltering ? cfg.getParameter<double>("phiPreFiltering") : 0) {
41  edm::ParameterSet comparitorPSet = cfg.getParameter<edm::ParameterSet>("SeedComparitorPSet");
42  std::string comparitorName = comparitorPSet.getParameter<std::string>("ComponentName");
43  if (comparitorName != "none") {
44  theComparitor = SeedComparitorFactory::get()->create(comparitorName, comparitorPSet, iC);
45  }
46  if (useMScat) {
48  }
49 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::unique_ptr< SeedComparitor > theComparitor
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > msmakerToken_
#define get
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > fieldToken_

◆ ~PixelTripletHLTGenerator()

PixelTripletHLTGenerator::~PixelTripletHLTGenerator ( )
override

Definition at line 51 of file PixelTripletHLTGenerator.cc.

51 {}

Member Function Documentation

◆ fillDescriptions()

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

Definition at line 53 of file PixelTripletHLTGenerator.cc.

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

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

◆ fillDescriptionsLabel()

static const char* PixelTripletHLTGenerator::fillDescriptionsLabel ( )
inlinestatic

Definition at line 34 of file PixelTripletHLTGenerator.h.

34 { return "pixelTripletHLT"; }

◆ hitTriplets() [1/4]

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 68 of file PixelTripletHLTGenerator.cc.

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

Referenced by hitTriplets().

73  {
74  auto const& doublets = thePairGenerator->doublets(region, ev, es, pairLayers);
75  if (not doublets or doublets->empty())
76  return;
77 
79  hitTriplets(region, result, ev, es, *doublets, thirdLayers, nullptr, *theLayerCache);
80 }
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
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

◆ hitTriplets() [2/4]

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 82 of file PixelTripletHLTGenerator.cc.

References HLT_2024v14_cff::doublets, makeMEIFBenchmarkPlots::ev, hitTriplets(), nano_mu_digi_cff::region, mps_fire::result, findQualityFiles::size, and theComparitor.

89  {
90  if (theComparitor)
91  theComparitor->init(ev, es);
92 
93  int size = thirdLayers.size();
94  const RecHitsSortedInPhi* thirdHitMap[size];
95  vector<const DetLayer*> thirdLayerDetLayer(size, nullptr);
96  for (int il = 0; il < size; ++il) {
97  thirdHitMap[il] = &layerCache(thirdLayers[il], region);
98  thirdLayerDetLayer[il] = thirdLayers[il].detLayer();
99  }
100  hitTriplets(region, result, es, doublets, thirdHitMap, thirdLayerDetLayer, size, tripletLastLayerIndex);
101 }
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

◆ hitTriplets() [3/4]

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 103 of file PixelTripletHLTGenerator.cc.

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

109  {
110  hitTriplets(region, result, es, doublets, thirdHitMap, thirdLayerDetLayer, nThirdLayers, nullptr);
111 }
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 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 113 of file PixelTripletHLTGenerator.cc.

References angle(), checkPhiInRange(), ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), pfMETCorrectionType0_cfi::correction, PixelRecoUtilities::curvature(), declareDynArray, HLT_2024v14_cff::doublets, dphi, relativeConstraints::empty, extraHitRPhitolerance, extraHitRZtolerance, f, fieldToken_, nano_mu_digi_cff::float, Geom::ftwoPi(), edm::EventSetup::getData(), hfClusterShapes_cfi::hits, mps_fire::i, ThirdHitRZPredictionBase::initLayer(), ThirdHitRZPredictionBase::initTolerance(), HitDoublets::inner, PixelRecoRange< T >::intersection(), nano_mu_digi_cff::layer, mps_splice::line, LogDebug, M_PI, SiStripPI::max, SiStripPI::min, msmakerToken_, normalizedPhi(), nSigmaPhi, nSigmaRZ, HitDoublets::outer, proxim(), CosmicsPD_Skims::radius, nano_mu_digi_cff::region, mps_fire::result, edm::swap(), theComparitor, HitTripletGeneratorFromPairAndLayers::theMaxElement, useBend, useFixedPreFiltering, useMScat, findQualityFiles::v, and protons_cff::xi.

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

Member Data Documentation

◆ dphi

const float PixelTripletHLTGenerator::dphi
private

Definition at line 77 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

◆ extraHitRPhitolerance

const float PixelTripletHLTGenerator::extraHitRPhitolerance
private

Definition at line 74 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

◆ extraHitRZtolerance

const float PixelTripletHLTGenerator::extraHitRZtolerance
private

Definition at line 73 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

◆ fieldToken_

const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> PixelTripletHLTGenerator::fieldToken_
private

Definition at line 70 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

◆ msmakerToken_

Definition at line 71 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets(), and PixelTripletHLTGenerator().

◆ theComparitor

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

Definition at line 78 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets(), and PixelTripletHLTGenerator().

◆ useBend

const bool PixelTripletHLTGenerator::useBend
private

Definition at line 76 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

◆ useFixedPreFiltering

const bool PixelTripletHLTGenerator::useFixedPreFiltering
private

Definition at line 72 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

◆ useMScat

const bool PixelTripletHLTGenerator::useMScat
private

Definition at line 75 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets(), and PixelTripletHLTGenerator().