CMS 3D CMS Logo

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