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 using namespace ctfseeding;
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 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  const std::vector<SeedingLayer> &layers,
59  LayerCacheType *layerCache)
60 {
61  thePairGenerator = pairs.clone();
62  theLayers = layers;
63  theLayerCache = layerCache;
64 }
65 
66 
67 namespace {
68  inline
69  bool intersect(Range &range, const Range &second)
70  {
71  if (range.first > second.max() || range.second < second.min())
72  return false;
73  if (range.first < second.min())
74  range.first = second.min();
75  if (range.second > second.max())
76  range.second = second.max();
77  return range.first < range.second;
78  }
79 }
80 
83  const edm::Event & ev,
84  const edm::EventSetup& es)
85 {
87  es.get<TrackerDigiGeometryRecord>().get(tracker);
88 
89  //Retrieve tracker topology from geometry
91  es.get<IdealGeometryRecord>().get(tTopoHand);
92  const TrackerTopology *tTopo=tTopoHand.product();
93 
94  auto const & doublets = thePairGenerator->doublets(region,ev,es);
95 
96  if (doublets.empty()) return;
97 
98  auto outSeq = doublets.detLayer(HitDoublets::outer)->seqNum();
99 
100 
101  int size = theLayers.size();
102 
103 
104  using NodeInfo = KDTreeNodeInfo<unsigned int>;
105  std::vector<NodeInfo > layerTree; // re-used throughout
107 
108  float rzError[size]; //save maximum errors
109  float maxphi = Geom::ftwoPi(), minphi = -maxphi; //increase to cater for any range
110 
111  LayerRZPredictions mapPred[size];
112 
113  const RecHitsSortedInPhi * thirdHitMap[size];
114 
115  for(int il = 0; il < size; il++) {
116  thirdHitMap[il] = &(*theLayerCache)(&theLayers[il], region, ev, es);
117  auto const & hits = *thirdHitMap[il];
118 
119  const DetLayer *layer = theLayers[il].detLayer();
120  LayerRZPredictions &predRZ = mapPred[il];
121  predRZ.line.initLayer(layer);
122  predRZ.helix1.initLayer(layer);
123  predRZ.helix2.initLayer(layer);
124  predRZ.line.initTolerance(extraHitRZtolerance);
125  predRZ.helix1.initTolerance(extraHitRZtolerance);
126  predRZ.helix2.initTolerance(extraHitRZtolerance);
127  predRZ.rzPositionFixup = MatchedHitRZCorrectionFromBending(layer,tTopo);
128 
129  layerTree.clear();
130  float minv=999999.0; float maxv = -999999.0; // Initialise to extreme values in case no hits
131  float maxErr=0.0f;
132  for (unsigned int i=0; i!=hits.size(); ++i) {
133  auto angle = hits.phi(i);
134  auto v = hits.gv(i);
135  //use (phi,r) for endcaps rather than (phi,z)
136  minv = std::min(minv,v); maxv = std::max(maxv,v);
137  float myerr = hits.dv[i];
138  maxErr = std::max(maxErr,myerr);
139  layerTree.emplace_back(i, angle, v); // save it
140  if (angle < 0) // wrap all points in phi
141  { layerTree.emplace_back(i, angle+Geom::ftwoPi(), v);}
142  else
143  { layerTree.emplace_back(i, angle-Geom::ftwoPi(), v);}
144  }
145  KDTreeBox phiZ(minphi, maxphi, minv-0.01f, maxv+0.01f); // declare our bounds
146  //add fudge factors in case only one hit and also for floating-point inaccuracy
147  hitTree[il].build(layerTree, phiZ); // make KDtree
148  rzError[il] = maxErr; //save error
149  }
150 
151  double curv = PixelRecoUtilities::curvature(1. / region.ptMin(), es);
152 
153  for (std::size_t ip =0; ip!=doublets.size(); ip++) {
154  auto xi = doublets.x(ip,HitDoublets::inner);
155  auto yi = doublets.y(ip,HitDoublets::inner);
156  auto zi = doublets.z(ip,HitDoublets::inner);
157  // auto rvi = doublets.rv(ip,HitDoublets::inner);
158  auto xo = doublets.x(ip,HitDoublets::outer);
159  auto yo = doublets.y(ip,HitDoublets::outer);
160  auto zo = doublets.z(ip,HitDoublets::outer);
161  // auto rvo = doublets.rv(ip,HitDoublets::outer);
162  GlobalPoint gp1(xi,yi,zi);
163  GlobalPoint gp2(xo,yo,zo);
164 
165  PixelRecoLineRZ line(gp1, gp2);
166  PixelRecoPointRZ point2(gp2.perp(), zo);
167  ThirdHitPredictionFromCircle predictionRPhi(gp1, gp2, extraHitRPhitolerance);
168 
169  Range generalCurvature = predictionRPhi.curvature(region.originRBound());
170  if (!intersect(generalCurvature, Range(-curv, curv))) continue;
171 
172  for(int il = 0; il < size; il++) {
173  if (hitTree[il].empty()) continue; // Don't bother if no hits
174  const DetLayer *layer = theLayers[il].detLayer();
175  bool barrelLayer = layer->isBarrel();
176 
177  Range curvature = generalCurvature;
178  ThirdHitCorrection correction(es, region.ptMin(), layer, line, point2, outSeq, useMScat);
179 
180  LayerRZPredictions &predRZ = mapPred[il];
181  predRZ.line.initPropagator(&line);
182 
183  Range rzRange;
184  if (useBend) {
185  // For the barrel region:
186  // swiping the helix passing through the two points across from
187  // negative to positive bending, can give us a sort of U-shaped
188  // projection onto the phi-z (barrel) or r-z plane (forward)
189  // so we checking minimum/maximum of all three possible extrema
190  //
191  // For the endcap region:
192  // Checking minimum/maximum radius of the helix projection
193  // onto an endcap plane, here we have to guard against
194  // looping tracks, when phi(delta z) gets out of control.
195  // HelixRZ::rAtZ should not follow looping tracks, but clamp
196  // to the minimum reachable r with the next-best lower |curvature|.
197  // So same procedure as for the barrel region can be applied.
198  //
199  // In order to avoid looking for potential looping tracks at all
200  // we also clamp the allowed curvature range for this layer,
201  // and potentially fail the layer entirely
202 
203  if (!barrelLayer) {
204  Range z3s = predRZ.line.detRange();
205  double z3 = z3s.first < 0 ? std::max(z3s.first, z3s.second)
206  : std::min(z3s.first, z3s.second);
207  double maxCurvature = HelixRZ::maxCurvature(&predictionRPhi,
208  gp1.z(), gp2.z(), z3);
209  if (!intersect(curvature, Range(-maxCurvature, maxCurvature)))
210  continue;
211  }
212 
213  HelixRZ helix1(&predictionRPhi, gp1.z(), gp2.z(), curvature.first);
214  HelixRZ helix2(&predictionRPhi, gp1.z(), gp2.z(), curvature.second);
215 
216  predRZ.helix1.initPropagator(&helix1);
217  predRZ.helix2.initPropagator(&helix2);
218 
219  Range rzRanges[2] = { predRZ.helix1(), predRZ.helix2() };
220  predRZ.helix1.initPropagator(nullptr);
221  predRZ.helix2.initPropagator(nullptr);
222 
223  rzRange.first = std::min(rzRanges[0].first, rzRanges[1].first);
224  rzRange.second = std::max(rzRanges[0].second, rzRanges[1].second);
225 
226  // if the allowed curvatures include a straight line,
227  // this can give us another extremum for allowed r/z
228  if (curvature.first * curvature.second < 0.0) {
229  Range rzLineRange = predRZ.line();
230  rzRange.first = std::min(rzRange.first, rzLineRange.first);
231  rzRange.second = std::max(rzRange.second, rzLineRange.second);
232  }
233  } else {
234  rzRange = predRZ.line();
235  }
236 
237  if (rzRange.first >= rzRange.second)
238  continue;
239 
240  correction.correctRZRange(rzRange);
241 
242  Range phiRange;
243  if (useFixedPreFiltering) {
244  float phi0 = doublets.phi(ip,HitDoublets::outer);
245  phiRange = Range(phi0 - dphi, phi0 + dphi);
246  } else {
247  Range radius;
248 
249  if (barrelLayer) {
250  radius = predRZ.line.detRange();
251  if (!intersect(rzRange, predRZ.line.detSize()))
252  continue;
253  } else {
254  radius = rzRange;
255  if (!intersect(radius, predRZ.line.detSize()))
256  continue;
257  }
258 
259  Range rPhi1 = predictionRPhi(curvature, radius.first);
260  Range rPhi2 = predictionRPhi(curvature, radius.second);
261  correction.correctRPhiRange(rPhi1);
262  correction.correctRPhiRange(rPhi2);
263  rPhi1.first /= radius.first;
264  rPhi1.second /= radius.first;
265  rPhi2.first /= radius.second;
266  rPhi2.second /= radius.second;
267  phiRange = mergePhiRanges(rPhi1, rPhi2);
268  }
269 
270  layerTree.clear(); // Now recover hits in bounding box...
271  float prmin=phiRange.min(), prmax=phiRange.max(); //get contiguous range
272  if ((prmax-prmin) > Geom::ftwoPi())
273  { prmax=Geom::fpi(); prmin = -Geom::fpi();}
274  else
275  { while (prmax>maxphi) { prmin -= Geom::ftwoPi(); prmax -= Geom::ftwoPi();}
276  while (prmin<minphi) { prmin += Geom::ftwoPi(); prmax += Geom::ftwoPi();}
277  // This needs range -twoPi to +twoPi to work
278  }
279  if (barrelLayer) {
280  Range regMax = predRZ.line.detRange();
281  Range regMin = predRZ.line(regMax.min());
282  regMax = predRZ.line(regMax.max());
283  correction.correctRZRange(regMin);
284  correction.correctRZRange(regMax);
285  if (regMax.min() < regMin.min()) { swap(regMax, regMin);}
286  KDTreeBox phiZ(prmin, prmax,
287  regMin.min()-fnSigmaRZ*rzError[il],
288  regMax.max()+fnSigmaRZ*rzError[il]);
289  hitTree[il].search(phiZ, layerTree);
290  }
291  else {
292  KDTreeBox phiZ(prmin, prmax,
293  rzRange.min()-fnSigmaRZ*rzError[il],
294  rzRange.max()+fnSigmaRZ*rzError[il]);
295  hitTree[il].search(phiZ, layerTree);
296  }
297 
298  MatchedHitRZCorrectionFromBending l2rzFixup(doublets.hit(ip,HitDoublets::outer)->det()->geographicalId(), tTopo);
299  MatchedHitRZCorrectionFromBending l3rzFixup = predRZ.rzPositionFixup;
300 
301  thirdHitMap[il] = &(*theLayerCache)(&theLayers[il], region, ev, es);
302  auto const & hits = *thirdHitMap[il];
303  for (auto const & ih : layerTree) {
304  auto KDdata = ih.data;
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
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
#define min(a, b)
Definition: mlp_lapack.h:161
T r() const
Radius, same as mag()
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)
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]
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:94
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
std::vector< ctfseeding::SeedingLayer > theLayers
T const * product() const
Definition: ESHandle.h:62
float ptMin() const
minimal pt of interest
PixelTripletLargeTipGenerator(const edm::ParameterSet &cfg)
virtual void init(const HitPairGenerator &pairs, const std::vector< ctfseeding::SeedingLayer > &layers, LayerCacheType *layerCache)
virtual unsigned int size() const
T x() const
Definition: PV3DBase.h:62
DetLayer const * detLayer(layer l) const
tuple size
Write out results.
#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