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
106  std::vector<unsigned int> foundNodes; // re-used throughout
107  foundNodes.reserve(100);
109 
110  float rzError[size]; //save maximum errors
111  float maxphi = Geom::ftwoPi(), minphi = -maxphi; //increase to cater for any range
112 
113  LayerRZPredictions mapPred[size];
114 
115  const RecHitsSortedInPhi * thirdHitMap[size];
116 
117  for(int il = 0; il < size; il++) {
118  thirdHitMap[il] = &(*theLayerCache)(&theLayers[il], region, ev, es);
119  auto const & hits = *thirdHitMap[il];
120 
121  const DetLayer *layer = theLayers[il].detLayer();
122  LayerRZPredictions &predRZ = mapPred[il];
123  predRZ.line.initLayer(layer);
124  predRZ.helix1.initLayer(layer);
125  predRZ.helix2.initLayer(layer);
126  predRZ.line.initTolerance(extraHitRZtolerance);
127  predRZ.helix1.initTolerance(extraHitRZtolerance);
128  predRZ.helix2.initTolerance(extraHitRZtolerance);
129  predRZ.rzPositionFixup = MatchedHitRZCorrectionFromBending(layer,tTopo);
130 
131  layerTree.clear();
132  float minv=999999.0; float maxv = -999999.0; // Initialise to extreme values in case no hits
133  float maxErr=0.0f;
134  for (unsigned int i=0; i!=hits.size(); ++i) {
135  auto angle = hits.phi(i);
136  auto v = hits.gv(i);
137  //use (phi,r) for endcaps rather than (phi,z)
138  minv = std::min(minv,v); maxv = std::max(maxv,v);
139  float myerr = hits.dv[i];
140  maxErr = std::max(maxErr,myerr);
141  layerTree.emplace_back(i, angle, v); // save it
142  if (angle < 0) // wrap all points in phi
143  { layerTree.emplace_back(i, angle+Geom::ftwoPi(), v);}
144  else
145  { layerTree.emplace_back(i, angle-Geom::ftwoPi(), v);}
146  }
147  KDTreeBox phiZ(minphi, maxphi, minv-0.01f, maxv+0.01f); // declare our bounds
148  //add fudge factors in case only one hit and also for floating-point inaccuracy
149  hitTree[il].build(layerTree, phiZ); // make KDtree
150  rzError[il] = maxErr; //save error
151  }
152 
153  double curv = PixelRecoUtilities::curvature(1. / region.ptMin(), es);
154 
155  for (std::size_t ip =0; ip!=doublets.size(); ip++) {
156  auto xi = doublets.x(ip,HitDoublets::inner);
157  auto yi = doublets.y(ip,HitDoublets::inner);
158  auto zi = doublets.z(ip,HitDoublets::inner);
159  // auto rvi = doublets.rv(ip,HitDoublets::inner);
160  auto xo = doublets.x(ip,HitDoublets::outer);
161  auto yo = doublets.y(ip,HitDoublets::outer);
162  auto zo = doublets.z(ip,HitDoublets::outer);
163  // auto rvo = doublets.rv(ip,HitDoublets::outer);
164  GlobalPoint gp1(xi,yi,zi);
165  GlobalPoint gp2(xo,yo,zo);
166 
167  PixelRecoLineRZ line(gp1, gp2);
168  PixelRecoPointRZ point2(gp2.perp(), zo);
169  ThirdHitPredictionFromCircle predictionRPhi(gp1, gp2, extraHitRPhitolerance);
170 
171  Range generalCurvature = predictionRPhi.curvature(region.originRBound());
172  if (!intersect(generalCurvature, Range(-curv, curv))) continue;
173 
174  for(int il = 0; il < size; il++) {
175  if (hitTree[il].empty()) continue; // Don't bother if no hits
176  const DetLayer *layer = theLayers[il].detLayer();
177  bool barrelLayer = layer->isBarrel();
178 
179  Range curvature = generalCurvature;
180  ThirdHitCorrection correction(es, region.ptMin(), layer, line, point2, outSeq, useMScat);
181 
182  LayerRZPredictions &predRZ = mapPred[il];
183  predRZ.line.initPropagator(&line);
184 
185  Range rzRange;
186  if (useBend) {
187  // For the barrel region:
188  // swiping the helix passing through the two points across from
189  // negative to positive bending, can give us a sort of U-shaped
190  // projection onto the phi-z (barrel) or r-z plane (forward)
191  // so we checking minimum/maximum of all three possible extrema
192  //
193  // For the endcap region:
194  // Checking minimum/maximum radius of the helix projection
195  // onto an endcap plane, here we have to guard against
196  // looping tracks, when phi(delta z) gets out of control.
197  // HelixRZ::rAtZ should not follow looping tracks, but clamp
198  // to the minimum reachable r with the next-best lower |curvature|.
199  // So same procedure as for the barrel region can be applied.
200  //
201  // In order to avoid looking for potential looping tracks at all
202  // we also clamp the allowed curvature range for this layer,
203  // and potentially fail the layer entirely
204 
205  if (!barrelLayer) {
206  Range z3s = predRZ.line.detRange();
207  double z3 = z3s.first < 0 ? std::max(z3s.first, z3s.second)
208  : std::min(z3s.first, z3s.second);
209  double maxCurvature = HelixRZ::maxCurvature(&predictionRPhi,
210  gp1.z(), gp2.z(), z3);
211  if (!intersect(curvature, Range(-maxCurvature, maxCurvature)))
212  continue;
213  }
214 
215  HelixRZ helix1(&predictionRPhi, gp1.z(), gp2.z(), curvature.first);
216  HelixRZ helix2(&predictionRPhi, gp1.z(), gp2.z(), curvature.second);
217 
218  predRZ.helix1.initPropagator(&helix1);
219  predRZ.helix2.initPropagator(&helix2);
220 
221  Range rzRanges[2] = { predRZ.helix1(), predRZ.helix2() };
222  predRZ.helix1.initPropagator(nullptr);
223  predRZ.helix2.initPropagator(nullptr);
224 
225  rzRange.first = std::min(rzRanges[0].first, rzRanges[1].first);
226  rzRange.second = std::max(rzRanges[0].second, rzRanges[1].second);
227 
228  // if the allowed curvatures include a straight line,
229  // this can give us another extremum for allowed r/z
230  if (curvature.first * curvature.second < 0.0) {
231  Range rzLineRange = predRZ.line();
232  rzRange.first = std::min(rzRange.first, rzLineRange.first);
233  rzRange.second = std::max(rzRange.second, rzLineRange.second);
234  }
235  } else {
236  rzRange = predRZ.line();
237  }
238 
239  if (rzRange.first >= rzRange.second)
240  continue;
241 
242  correction.correctRZRange(rzRange);
243 
244  Range phiRange;
245  if (useFixedPreFiltering) {
246  float phi0 = doublets.phi(ip,HitDoublets::outer);
247  phiRange = Range(phi0 - dphi, phi0 + dphi);
248  } else {
249  Range radius;
250 
251  if (barrelLayer) {
252  radius = predRZ.line.detRange();
253  if (!intersect(rzRange, predRZ.line.detSize()))
254  continue;
255  } else {
256  radius = rzRange;
257  if (!intersect(radius, predRZ.line.detSize()))
258  continue;
259  }
260 
261  Range rPhi1 = predictionRPhi(curvature, radius.first);
262  Range rPhi2 = predictionRPhi(curvature, radius.second);
263  correction.correctRPhiRange(rPhi1);
264  correction.correctRPhiRange(rPhi2);
265  rPhi1.first /= radius.first;
266  rPhi1.second /= radius.first;
267  rPhi2.first /= radius.second;
268  rPhi2.second /= radius.second;
269  phiRange = mergePhiRanges(rPhi1, rPhi2);
270  }
271 
272  foundNodes.clear(); // Now recover hits in bounding box...
273  float prmin=phiRange.min(), prmax=phiRange.max(); //get contiguous range
274  if ((prmax-prmin) > Geom::ftwoPi())
275  { prmax=Geom::fpi(); prmin = -Geom::fpi();}
276  else
277  { while (prmax>maxphi) { prmin -= Geom::ftwoPi(); prmax -= Geom::ftwoPi();}
278  while (prmin<minphi) { prmin += Geom::ftwoPi(); prmax += Geom::ftwoPi();}
279  // This needs range -twoPi to +twoPi to work
280  }
281  if (barrelLayer) {
282  Range regMax = predRZ.line.detRange();
283  Range regMin = predRZ.line(regMax.min());
284  regMax = predRZ.line(regMax.max());
285  correction.correctRZRange(regMin);
286  correction.correctRZRange(regMax);
287  if (regMax.min() < regMin.min()) { swap(regMax, regMin);}
288  KDTreeBox phiZ(prmin, prmax,
289  regMin.min()-fnSigmaRZ*rzError[il],
290  regMax.max()+fnSigmaRZ*rzError[il]);
291  hitTree[il].search(phiZ, foundNodes);
292  }
293  else {
294  KDTreeBox phiZ(prmin, prmax,
295  rzRange.min()-fnSigmaRZ*rzError[il],
296  rzRange.max()+fnSigmaRZ*rzError[il]);
297  hitTree[il].search(phiZ, foundNodes);
298  }
299 
300  MatchedHitRZCorrectionFromBending l2rzFixup(doublets.hit(ip,HitDoublets::outer)->det()->geographicalId(), tTopo);
301  MatchedHitRZCorrectionFromBending l3rzFixup = predRZ.rzPositionFixup;
302 
303  thirdHitMap[il] = &(*theLayerCache)(&theLayers[il], region, ev, es);
304  auto const & hits = *thirdHitMap[il];
305  for (auto KDdata : foundNodes) {
306  GlobalPoint p3 = hits.gp(KDdata);
307  double p3_r = p3.perp();
308  double p3_z = p3.z();
309  float p3_phi = hits.phi(KDdata);
310 
311  Range rangeRPhi = predictionRPhi(curvature, p3_r);
312  correction.correctRPhiRange(rangeRPhi);
313 
314  float ir = 1.f/p3_r;
315  float phiErr = nSigmaPhi * hits.drphi[KDdata]*ir;
316  if (!checkPhiInRange(p3_phi, rangeRPhi.first*ir-phiErr, rangeRPhi.second*ir+phiErr))
317  continue;
318 
319  Basic2DVector<double> thc(p3.x(), p3.y());
320 
321  auto curv_ = predictionRPhi.curvature(thc);
322  double p2_r = point2.r(); double p2_z = point2.z(); // they will be modified!
323 
324  l2rzFixup(predictionRPhi, curv_, *doublets.hit(ip,HitDoublets::outer), p2_r, p2_z, tTopo);
325  l3rzFixup(predictionRPhi, curv_, *hits.theHits[KDdata].hit(), p3_r, p3_z, tTopo);
326 
327  Range rangeRZ;
328  if (useBend) {
329  HelixRZ updatedHelix(&predictionRPhi, gp1.z(), p2_z, curv_);
330  rangeRZ = predRZ.helix1(barrelLayer ? p3_r : p3_z, updatedHelix);
331  } else {
332  float tIP = predictionRPhi.transverseIP(thc);
333  PixelRecoPointRZ updatedPoint2(p2_r, p2_z);
334  PixelRecoLineRZ updatedLine(line.origin(), point2, tIP);
335  rangeRZ = predRZ.line(barrelLayer ? p3_r : p3_z, line);
336  }
337  correction.correctRZRange(rangeRZ);
338 
339  double err = nSigmaRZ * hits.dv[KDdata];
340 
341  rangeRZ.first -= err, rangeRZ.second += err;
342 
343  if (!rangeRZ.inside(barrelLayer ? p3_z : p3_r)) continue;
344 
345  if (theMaxElement!=0 && result.size() >= theMaxElement) {
346  result.clear();
347  edm::LogError("TooManyTriplets")<<" number of triples exceed maximum. no triplets produced.";
348  return;
349  }
350  result.emplace_back( doublets.hit(ip,HitDoublets::inner), doublets.hit(ip,HitDoublets::outer), hits.theHits[KDdata].hit());
351  }
352  }
353  }
354  // std::cout << "found triplets " << result.size() << std::endl;
355 }
356 
357 bool PixelTripletLargeTipGenerator::checkPhiInRange(float phi, float phi1, float phi2) const
358 { while (phi > phi2) phi -= 2. * M_PI;
359  while (phi < phi1) phi += 2. * M_PI;
360  return phi <= phi2;
361 }
362 
363 std::pair<float, float>
364 PixelTripletLargeTipGenerator::mergePhiRanges(const std::pair<float, float> &r1,
365  const std::pair<float, float> &r2) const
366 { float r2Min = r2.first;
367  float r2Max = r2.second;
368  while (r1.first - r2Min > +M_PI) r2Min += 2. * M_PI, r2Max += 2. * M_PI;
369  while (r1.first - r2Min < -M_PI) r2Min -= 2. * M_PI, r2Max -= 2. * M_PI;
370  return std::make_pair(min(r1.first, r2Min), max(r1.second, r2Max));
371 }
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()
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: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
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