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