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  theComparitor = (comparitorName == "none") ?
48  0 : SeedComparitorFactory::get()->create( comparitorName, comparitorPSet);
49 }
50 
52  delete thePairGenerator;
53  delete theComparitor;
54 }
55 
57  const std::vector<SeedingLayer> & layers,
58  LayerCacheType* layerCache)
59 {
60  thePairGenerator = pairs.clone();
61  theLayers = layers;
62  theLayerCache = layerCache;
63 }
64 
67  const edm::Event & ev,
68  const edm::EventSetup& es)
69 {
70 
72 
73  auto const & doublets = thePairGenerator->doublets(region,ev,es);
74 
75  if (doublets.empty()) return;
76 
77  auto outSeq = doublets.detLayer(HitDoublets::outer)->seqNum();
78 
79 
80  // std::cout << "pairs " << doublets.size() << std::endl;
81 
82  float regOffset = region.origin().perp(); //try to take account of non-centrality (?)
83  int size = theLayers.size();
84 
86 
87  const RecHitsSortedInPhi * thirdHitMap[size];
89 
90  using NodeInfo = KDTreeNodeInfo<unsigned int>;
91  std::vector<NodeInfo > layerTree; // re-used throughout
92 
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)(&theLayers[il], region, ev, es);
100  auto const & hits = *thirdHitMap[il];
101  ThirdHitRZPrediction<PixelRecoLineRZ> & pred = preds[il];
102  pred.initLayer(theLayers[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 " << theLayers[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  ThirdHitPredictionFromInvParabola predictionRPhitmp(xi,yi,xo,yo,imppartmp,curv,extraHitRPhitolerance);
149 
150  // printf("++Constr %f %f %f %f %f %f %f\n",xi,yi,xo,yo,imppartmp,curv,extraHitRPhitolerance);
151 
152  // std::cout << ip << ": " << point1.r() << ","<< point1.z() << " "
153  // << point2.r() << ","<< point2.z() <<std::endl;
154 
155  for (int il=0; il!=size; ++il) {
156  if (hitTree[il].empty()) continue; // Don't bother if no hits
157 
158  auto const & hits = *thirdHitMap[il];
159 
160  const DetLayer * layer = theLayers[il].detLayer();
161  auto barrelLayer = layer->isBarrel();
162 
163  ThirdHitCorrection correction(es, region.ptMin(), layer, line, point2, outSeq, useMScat, useBend);
164 
165  ThirdHitRZPrediction<PixelRecoLineRZ> & predictionRZ = preds[il];
166 
167  predictionRZ.initPropagator(&line);
168  Range rzRange = predictionRZ();
169  correction.correctRZRange(rzRange);
170 
171  Range phiRange;
172  if (useFixedPreFiltering) {
173  float phi0 = doublets.phi(ip,HitDoublets::outer);
174  phiRange = Range(phi0-dphi,phi0+dphi);
175  }
176  else {
177  Range radius;
178  if (barrelLayer) {
179  radius = predictionRZ.detRange();
180  } else {
181  radius = Range(max(rzRange.min(), predictionRZ.detSize().min()),
182  min(rzRange.max(), predictionRZ.detSize().max()) );
183  }
184  if (radius.empty()) continue;
185 
186  // std::cout << "++R " << radius.min() << " " << radius.max() << std::endl;
187 
188  Range rPhi1m = predictionRPhitmp(radius.max(), -1);
189  Range rPhi1p = predictionRPhitmp(radius.max(), 1);
190  Range rPhi2m = predictionRPhitmp(radius.min(), -1);
191  Range rPhi2p = predictionRPhitmp(radius.min(), 1);
192  Range rPhi1 = rPhi1m.sum(rPhi1p);
193  Range rPhi2 = rPhi2m.sum(rPhi2p);
194  correction.correctRPhiRange(rPhi1);
195  correction.correctRPhiRange(rPhi2);
196  rPhi1.first /= radius.max();
197  rPhi1.second /= radius.max();
198  rPhi2.first /= radius.min();
199  rPhi2.second /= radius.min();
200  phiRange = mergePhiRanges(rPhi1,rPhi2);
201  }
202 
203  constexpr float nSigmaRZ = std::sqrt(12.f); // ...and continue as before
204  constexpr float nSigmaPhi = 3.f;
205 
206  layerTree.clear(); // Now recover hits in bounding box...
207  float prmin=phiRange.min(), prmax=phiRange.max();
208  if ((prmax-prmin) > Geom::ftwoPi())
209  { prmax=Geom::fpi(); prmin = -Geom::fpi();}
210  else
211  { while (prmax>maxphi) { prmin -= Geom::ftwoPi(); prmax -= Geom::ftwoPi();}
212  while (prmin<minphi) { prmin += Geom::ftwoPi(); prmax += Geom::ftwoPi();}
213  // This needs range -twoPi to +twoPi to work
214  }
215  if (barrelLayer)
216  {
217  Range regMax = predictionRZ.detRange();
218  Range regMin = predictionRZ(regMax.min()-regOffset);
219  regMax = predictionRZ(regMax.max()+regOffset);
220  correction.correctRZRange(regMin);
221  correction.correctRZRange(regMax);
222  if (regMax.min() < regMin.min()) { swap(regMax, regMin);}
223  KDTreeBox phiZ(prmin, prmax, regMin.min()-nSigmaRZ*rzError[il], regMax.max()+nSigmaRZ*rzError[il]);
224  hitTree[il].search(phiZ, layerTree);
225  }
226  else
227  {
228  KDTreeBox phiZ(prmin, prmax,
229  rzRange.min()-regOffset-nSigmaRZ*rzError[il],
230  rzRange.max()+regOffset+nSigmaRZ*rzError[il]);
231  hitTree[il].search(phiZ, layerTree);
232  }
233 
234  // std::cout << ip << ": " << theLayers[il].detLayer()->seqNum() << " " << layerTree.size() << " " << prmin << " " << prmax << std::endl;
235 
236 
237  // int kk=0;
238  for (auto const & ih : layerTree) {
239 
240  if (theMaxElement!=0 && result.size() >= theMaxElement){
241  result.clear();
242  edm::LogError("TooManyTriplets")<<" number of triples exceeds maximum. no triplets produced.";
243  return;
244  }
245 
246  auto KDdata = ih.data;
247  float p3_u = hits.u[KDdata];
248  float p3_v = hits.v[KDdata];
249  float p3_phi = hits.lphi[KDdata];
250 
251  //if ((kk++)%100==0)
252  //std::cout << kk << ": " << p3_u << " " << p3_v << " " << p3_phi << std::endl;
253 
254 
255  Range allowed = predictionRZ(p3_u);
256  correction.correctRZRange(allowed);
257  float vErr = nSigmaRZ *hits.dv[KDdata];
258  Range hitRange(p3_v-vErr, p3_v+vErr);
259  Range crossingRange = allowed.intersection(hitRange);
260  if (crossingRange.empty()) continue;
261 
262  float ir = 1.f/hits.rv(KDdata);
263  float phiErr = nSigmaPhi * hits.drphi[KDdata]*ir;
264  for (int icharge=-1; icharge <=1; icharge+=2) {
265  Range rangeRPhi = predictionRPhi(hits.rv(KDdata), icharge);
266  correction.correctRPhiRange(rangeRPhi);
267  if (checkPhiInRange(p3_phi, rangeRPhi.first*ir-phiErr, rangeRPhi.second*ir+phiErr)) {
268  // insert here check with comparitor
269  OrderedHitTriplet hittriplet( doublets.hit(ip,HitDoublets::inner), doublets.hit(ip,HitDoublets::outer), hits.theHits[KDdata].hit());
270  if (!theComparitor || theComparitor->compatible(hittriplet,region) ) {
271  result.push_back( hittriplet );
272  } else {
273  LogDebug("RejectedTriplet") << "rejected triplet from comparitor ";
274  }
275  break;
276  }
277  }
278  }
279  }
280  }
281  // std::cout << "triplets " << result.size() << std::endl;
282 }
283 
284 bool PixelTripletHLTGenerator::checkPhiInRange(float phi, float phi1, float phi2) const
285 {
286  while (phi > phi2) phi -= Geom::ftwoPi();
287  while (phi < phi1) phi += Geom::ftwoPi();
288  return ( (phi1 <= phi) && (phi <= phi2) );
289 }
290 
291 std::pair<float,float> PixelTripletHLTGenerator::mergePhiRanges(const std::pair<float,float>& r1,
292  const std::pair<float,float>& r2) const
293 { float r2_min=r2.first;
294  float r2_max=r2.second;
295  while (r1.first-r2_min > Geom::fpi()) { r2_min += Geom::ftwoPi(); r2_max += Geom::ftwoPi();}
296  while (r1.first-r2_min < -Geom::fpi()) { r2_min -= Geom::ftwoPi(); r2_max -= Geom::ftwoPi(); }
297 
298  return std::make_pair(min(r1.first,r2_min),max(r1.second,r2_max));
299 }
#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
std::vector< ctfseeding::SeedingLayer > theLayers
virtual HitPairGenerator * clone() const =0
void build(std::vector< KDTreeNodeInfo > &eltList, const KDTreeBox &region)
void initPropagator(const Propagator *propagator)
T perp() const
Definition: PV3DBase.h:72
GlobalPoint const & origin() const
T max() const
void initLayer(const DetLayer *layer)
PixelRecoRange< T > sum(const PixelRecoRange< T > &r) const
T y() const
Definition: PV3DBase.h:63
#define min(a, b)
Definition: mlp_lapack.h:161
virtual bool compatible(const SeedingHitSet &hits, const TrackingRegion &region) const =0
float fpi()
Definition: Pi.h:35
virtual void init(const edm::EventSetup &es)=0
virtual void init(const HitPairGenerator &pairs, const std::vector< ctfseeding::SeedingLayer > &layers, LayerCacheType *layerCache)
void search(const KDTreeBox &searchBox, std::vector< KDTreeNodeInfo > &resRecHitList)
std::pair< float, float > mergePhiRanges(const std::pair< float, float > &r1, const std::pair< float, float > &r2) const
int seqNum() const
Definition: DetLayer.h:39
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
tuple result
Definition: query.py:137
double f[11][100]
virtual HitDoublets doublets(const TrackingRegion &reg, const edm::Event &ev, const edm::EventSetup &es)
constexpr double nSigmaPhi
bool isBarrel() const
Definition: DetLayer.h:35
PixelRecoRange< float > Range
PixelTripletHLTGenerator(const edm::ParameterSet &cfg)
float ptMin() const
minimal pt of interest
PixelRecoRange< T > intersection(const PixelRecoRange< T > &r) const
bool checkPhiInRange(float phi, float phi1, float phi2) const
TransientTrackingRecHit::ConstRecHitPointer Hit
void initTolerance(float tolerance)
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
tuple size
Write out results.
T get(const Candidate &c)
Definition: component.h:56
#define constexpr
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