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