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