CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelTripletHLTGenerator.cc
Go to the documentation of this file.
3 
5 #include "ThirdHitRZPrediction.h"
8 
11 #include "ThirdHitCorrection.h"
14 #include <iostream>
15 
18 
20 //#include "RecoParticleFlow/PFProducer/interface/KDTreeLinkerAlgo.h"
21 //#include "RecoParticleFlow/PFProducer/interface/KDTreeLinkerTools.h"
22 #include "RecoPixelVertexing/PixelTriplets/plugins/KDTreeLinkerAlgo.h" //amend to point at your copy...
24 
25 #include<cstdio>
26 
29 
30 using namespace std;
31 using namespace ctfseeding;
32 
35  useFixedPreFiltering(cfg.getParameter<bool>("useFixedPreFiltering")),
36  extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
37  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
38  useMScat(cfg.getParameter<bool>("useMultScattering")),
39  useBend(cfg.getParameter<bool>("useBending"))
40 {
41  dphi = (useFixedPreFiltering) ? cfg.getParameter<double>("phiPreFiltering") : 0;
42 
43  edm::ParameterSet comparitorPSet =
44  cfg.getParameter<edm::ParameterSet>("SeedComparitorPSet");
45  std::string comparitorName = comparitorPSet.getParameter<std::string>("ComponentName");
46  if(comparitorName != "none") {
47  theComparitor.reset( SeedComparitorFactory::get()->create( comparitorName, comparitorPSet, iC) );
48  }
49 }
50 
52 
55  const edm::Event & ev,
56  const edm::EventSetup& es,
58  const std::vector<SeedingLayerSetsHits::SeedingLayer>& thirdLayers)
59 {
60 
61  if (theComparitor) theComparitor->init(ev, es);
62 
63  auto const & doublets = thePairGenerator->doublets(region,ev,es, pairLayers);
64 
65  if (doublets.empty()) return;
66 
67  auto outSeq = doublets.detLayer(HitDoublets::outer)->seqNum();
68 
69 
70  // std::cout << "pairs " << doublets.size() << std::endl;
71 
72  float regOffset = region.origin().perp(); //try to take account of non-centrality (?)
73  int size = thirdLayers.size();
74 
75  #ifdef __clang__
76  std::vector<ThirdHitRZPrediction<PixelRecoLineRZ>> preds(size);
77  #else
79  #endif
80 
81  const RecHitsSortedInPhi * thirdHitMap[size];
83 
84  using NodeInfo = KDTreeNodeInfo<unsigned int>;
85  std::vector<NodeInfo > layerTree; // re-used throughout
86  std::vector<unsigned int> foundNodes; // re-used thoughout
87  foundNodes.reserve(100);
88 
89  #ifdef __clang__
90  std::vector<KDTreeLinkerAlgo<unsigned int>> hitTree(size);
91  #else
93  #endif
94  float rzError[size]; //save maximum errors
95  float maxphi = Geom::ftwoPi(), minphi = -maxphi; // increase to cater for any range
96 
97  // fill the prediction vector
98  for (int il=0; il!=size; ++il) {
99  thirdHitMap[il] = &(*theLayerCache)(thirdLayers[il], region, ev, es);
100  auto const & hits = *thirdHitMap[il];
101  ThirdHitRZPrediction<PixelRecoLineRZ> & pred = preds[il];
102  pred.initLayer(thirdLayers[il].detLayer());
104 
105  layerTree.clear();
106  float minv=999999.0, maxv= -999999.0; // Initialise to extreme values in case no hits
107  float maxErr=0.0f;
108  for (unsigned int i=0; i!=hits.size(); ++i) {
109  auto angle = hits.phi(i);
110  auto v = hits.gv(i);
111  //use (phi,r) for endcaps rather than (phi,z)
112  minv = std::min(minv,v); maxv = std::max(maxv,v);
113  float myerr = hits.dv[i];
114  maxErr = std::max(maxErr,myerr);
115  layerTree.emplace_back(i, angle, v); // save it
116  if (angle < 0) // wrap all points in phi
117  { layerTree.emplace_back(i, angle+Geom::ftwoPi(), v);}
118  else
119  { layerTree.emplace_back(i, angle-Geom::ftwoPi(), v);}
120  }
121  KDTreeBox phiZ(minphi, maxphi, minv-0.01f, maxv+0.01f); // declare our bounds
122  //add fudge factors in case only one hit and also for floating-point inaccuracy
123  hitTree[il].build(layerTree, phiZ); // make KDtree
124  rzError[il] = maxErr; //save error
125  // std::cout << "layer " << thirdLayers[il].detLayer()->seqNum() << " " << layerTree.size() << std::endl;
126  }
127 
128  float imppar = region.originRBound();
129  float imppartmp = region.originRBound()+region.origin().perp();
130  float curv = PixelRecoUtilities::curvature(1.f/region.ptMin(), es);
131 
132  for (std::size_t ip =0; ip!=doublets.size(); ip++) {
133  auto xi = doublets.x(ip,HitDoublets::inner);
134  auto yi = doublets.y(ip,HitDoublets::inner);
135  auto zi = doublets.z(ip,HitDoublets::inner);
136  auto rvi = doublets.rv(ip,HitDoublets::inner);
137  auto xo = doublets.x(ip,HitDoublets::outer);
138  auto yo = doublets.y(ip,HitDoublets::outer);
139  auto zo = doublets.z(ip,HitDoublets::outer);
140  auto rvo = doublets.rv(ip,HitDoublets::outer);
141 
142  PixelRecoPointRZ point1(rvi, zi);
143  PixelRecoPointRZ point2(rvo, zo);
144  PixelRecoLineRZ line(point1, point2);
145  ThirdHitPredictionFromInvParabola predictionRPhi(xi-region.origin().x(),yi-region.origin().y(),
146  xo-region.origin().x(),yo-region.origin().y(),
147  imppar,curv,extraHitRPhitolerance);
148 
149  ThirdHitPredictionFromInvParabola predictionRPhitmp(xi,yi,xo,yo,imppartmp,curv,extraHitRPhitolerance);
150 
151  // printf("++Constr %f %f %f %f %f %f %f\n",xi,yi,xo,yo,imppartmp,curv,extraHitRPhitolerance);
152 
153  // std::cout << ip << ": " << point1.r() << ","<< point1.z() << " "
154  // << point2.r() << ","<< point2.z() <<std::endl;
155 
156  for (int il=0; il!=size; ++il) {
157  if (hitTree[il].empty()) continue; // Don't bother if no hits
158 
159  auto const & hits = *thirdHitMap[il];
160 
161  const DetLayer * layer = thirdLayers[il].detLayer();
162  auto barrelLayer = layer->isBarrel();
163 
164  ThirdHitCorrection correction(es, region.ptMin(), layer, line, point2, outSeq, useMScat, useBend);
165 
166  ThirdHitRZPrediction<PixelRecoLineRZ> & predictionRZ = preds[il];
167 
168  predictionRZ.initPropagator(&line);
169  Range rzRange = predictionRZ();
170  correction.correctRZRange(rzRange);
171 
172  Range phiRange;
173  if (useFixedPreFiltering) {
174  float phi0 = doublets.phi(ip,HitDoublets::outer);
175  phiRange = Range(phi0-dphi,phi0+dphi);
176  }
177  else {
178  Range radius;
179  if (barrelLayer) {
180  radius = predictionRZ.detRange();
181  } else {
182  radius = Range(max(rzRange.min(), predictionRZ.detSize().min()),
183  min(rzRange.max(), predictionRZ.detSize().max()) );
184  }
185  if (radius.empty()) continue;
186 
187  // std::cout << "++R " << radius.min() << " " << radius.max() << std::endl;
188 
189  /*
190  Range rPhi1m = predictionRPhitmp(radius.max(), -1);
191  Range rPhi1p = predictionRPhitmp(radius.max(), 1);
192  Range rPhi2m = predictionRPhitmp(radius.min(), -1);
193  Range rPhi2p = predictionRPhitmp(radius.min(), 1);
194  Range rPhi1 = rPhi1m.sum(rPhi1p);
195  Range rPhi2 = rPhi2m.sum(rPhi2p);
196 
197 
198  auto rPhi1N = predictionRPhitmp(radius.max());
199  auto rPhi2N = predictionRPhitmp(radius.min());
200 
201  std::cout << "VI "
202  << rPhi1N.first <<'/'<< rPhi1.first << ' '
203  << rPhi1N.second <<'/'<< rPhi1.second << ' '
204  << rPhi2N.first <<'/'<< rPhi2.first << ' '
205  << rPhi2N.second <<'/'<< rPhi2.second
206  << std::endl;
207 
208  */
209 
210  auto rPhi1 = predictionRPhitmp(radius.max());
211  auto rPhi2 = predictionRPhitmp(radius.min());
212 
213 
214  correction.correctRPhiRange(rPhi1);
215  correction.correctRPhiRange(rPhi2);
216  rPhi1.first /= radius.max();
217  rPhi1.second /= radius.max();
218  rPhi2.first /= radius.min();
219  rPhi2.second /= radius.min();
220  phiRange = mergePhiRanges(rPhi1,rPhi2);
221  }
222 
223  constexpr float nSigmaRZ = 3.46410161514f; // std::sqrt(12.f); // ...and continue as before
224  constexpr float nSigmaPhi = 3.f;
225 
226  foundNodes.clear(); // Now recover hits in bounding box...
227  float prmin=phiRange.min(), prmax=phiRange.max();
228  if ((prmax-prmin) > Geom::ftwoPi())
229  { prmax=Geom::fpi(); prmin = -Geom::fpi();}
230  else
231  { while (prmax>maxphi) { prmin -= Geom::ftwoPi(); prmax -= Geom::ftwoPi();}
232  while (prmin<minphi) { prmin += Geom::ftwoPi(); prmax += Geom::ftwoPi();}
233  // This needs range -twoPi to +twoPi to work
234  }
235  if (barrelLayer)
236  {
237  Range regMax = predictionRZ.detRange();
238  Range regMin = predictionRZ(regMax.min()-regOffset);
239  regMax = predictionRZ(regMax.max()+regOffset);
240  correction.correctRZRange(regMin);
241  correction.correctRZRange(regMax);
242  if (regMax.min() < regMin.min()) { swap(regMax, regMin);}
243  KDTreeBox phiZ(prmin, prmax, regMin.min()-nSigmaRZ*rzError[il], regMax.max()+nSigmaRZ*rzError[il]);
244  hitTree[il].search(phiZ, foundNodes);
245  }
246  else
247  {
248  KDTreeBox phiZ(prmin, prmax,
249  rzRange.min()-regOffset-nSigmaRZ*rzError[il],
250  rzRange.max()+regOffset+nSigmaRZ*rzError[il]);
251  hitTree[il].search(phiZ, foundNodes);
252  }
253 
254  // std::cout << ip << ": " << thirdLayers[il].detLayer()->seqNum() << " " << foundNodes.size() << " " << prmin << " " << prmax << std::endl;
255 
256 
257  // int kk=0;
258  for (auto KDdata : foundNodes) {
259 
260  if (theMaxElement!=0 && result.size() >= theMaxElement){
261  result.clear();
262  edm::LogError("TooManyTriplets")<<" number of triples exceeds maximum. no triplets produced.";
263  return;
264  }
265 
266  float p3_u = hits.u[KDdata];
267  float p3_v = hits.v[KDdata];
268  float p3_phi = hits.lphi[KDdata];
269 
270  //if ((kk++)%100==0)
271  //std::cout << kk << ": " << p3_u << " " << p3_v << " " << p3_phi << std::endl;
272 
273 
274  Range allowed = predictionRZ(p3_u);
275  correction.correctRZRange(allowed);
276  float vErr = nSigmaRZ *hits.dv[KDdata];
277  Range hitRange(p3_v-vErr, p3_v+vErr);
278  Range crossingRange = allowed.intersection(hitRange);
279  if (crossingRange.empty()) continue;
280 
281  float ir = 1.f/hits.rv(KDdata);
282  float phiErr = nSigmaPhi * hits.drphi[KDdata]*ir;
283  for (int icharge=-1; icharge <=1; icharge+=2) {
284  Range rangeRPhi = predictionRPhi(hits.rv(KDdata), icharge);
285  correction.correctRPhiRange(rangeRPhi);
286  if (checkPhiInRange(p3_phi, rangeRPhi.first*ir-phiErr, rangeRPhi.second*ir+phiErr)) {
287  // insert here check with comparitor
288  OrderedHitTriplet hittriplet( doublets.hit(ip,HitDoublets::inner), doublets.hit(ip,HitDoublets::outer), hits.theHits[KDdata].hit());
289  if (!theComparitor || theComparitor->compatible(hittriplet,region) ) {
290  result.push_back( hittriplet );
291  } else {
292  LogDebug("RejectedTriplet") << "rejected triplet from comparitor ";
293  }
294  break;
295  }
296  }
297  }
298  }
299  }
300  // std::cout << "triplets " << result.size() << std::endl;
301 }
302 
303 bool PixelTripletHLTGenerator::checkPhiInRange(float phi, float phi1, float phi2) const
304 {
305  while (phi > phi2) phi -= Geom::ftwoPi();
306  while (phi < phi1) phi += Geom::ftwoPi();
307  return ( (phi1 <= phi) && (phi <= phi2) );
308 }
309 
310 std::pair<float,float> PixelTripletHLTGenerator::mergePhiRanges(const std::pair<float,float>& r1,
311  const std::pair<float,float>& r2) const
312 { float r2_min=r2.first;
313  float r2_max=r2.second;
314  while (r1.first-r2_min > Geom::fpi()) { r2_min += Geom::ftwoPi(); r2_max += Geom::ftwoPi();}
315  while (r1.first-r2_min < -Geom::fpi()) { r2_min -= Geom::ftwoPi(); r2_max -= Geom::ftwoPi(); }
316 
317  return std::make_pair(min(r1.first,r2_min),max(r1.second,r2_max));
318 }
#define LogDebug(id)
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
void initPropagator(const Propagator *propagator)
T perp() const
Definition: PV3DBase.h:72
GlobalPoint const & origin() const
void initLayer(const DetLayer *layer)
T y() const
Definition: PV3DBase.h:63
bool ev
float fpi()
Definition: Pi.h:35
void search(const KDTreeBox &searchBox, std::vector< KDTreeNodeInfo > &resRecHitList)
#define constexpr
std::pair< float, float > mergePhiRanges(const std::pair< float, float > &r1, const std::pair< float, float > &r2) const
std::unique_ptr< SeedComparitor > theComparitor
PixelTripletHLTGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
T curvature(T InversePt, const edm::EventSetup &iSetup)
tuple result
Definition: query.py:137
double f[11][100]
T min(T a, T b)
Definition: MathUtil.h:58
SeedingHitSet::ConstRecHitPointer Hit
BaseTrackerRecHit const * Hit
bool isBarrel() const
Definition: DetLayer.h:32
PixelRecoRange< float > Range
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
std::unique_ptr< HitPairGeneratorFromLayerPair > thePairGenerator
float ptMin() const
minimal pt of interest
Geom::Phi< T > phi() const
PixelRecoRange< T > intersection(const PixelRecoRange< T > &r) const
bool checkPhiInRange(float phi, float phi1, float phi2) const
void initTolerance(float tolerance)
virtual unsigned int size() const
T x() const
Definition: PV3DBase.h:62
SurfaceDeformation * create(int type, const std::vector< double > &params)
tuple size
Write out results.
T get(const Candidate &c)
Definition: component.h:55
float ftwoPi()
Definition: Pi.h:36
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11