test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Attributes
PixelTripletLargeTipGenerator Class Reference

#include <PixelTripletLargeTipGenerator.h>

Inheritance diagram for PixelTripletLargeTipGenerator:
HitTripletGeneratorFromPairAndLayers

Public Member Functions

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
 
 PixelTripletLargeTipGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
 
virtual ~PixelTripletLargeTipGenerator ()
 
- Public Member Functions inherited from HitTripletGeneratorFromPairAndLayers
 HitTripletGeneratorFromPairAndLayers (unsigned int maxElement=0)
 
 HitTripletGeneratorFromPairAndLayers (const edm::ParameterSet &pset)
 
void init (std::unique_ptr< HitPairGeneratorFromLayerPair > &&pairs, LayerCacheType *layerCache)
 
const
HitPairGeneratorFromLayerPair
pairGenerator () const
 
virtual ~HitTripletGeneratorFromPairAndLayers ()
 

Private Types

typedef
CombinedHitTripletGenerator::LayerCacheType 
LayerCacheType
 

Private Attributes

const float dphi
 
const float extraHitRPhitolerance
 
const float extraHitRZtolerance
 
const bool useBend
 
const bool useFixedPreFiltering
 
const bool useMScat
 

Additional Inherited Members

- Public Types inherited from HitTripletGeneratorFromPairAndLayers
typedef LayerHitMapCache LayerCacheType
 
- Protected Attributes inherited from HitTripletGeneratorFromPairAndLayers
LayerCacheTypetheLayerCache
 
const unsigned int theMaxElement
 
std::unique_ptr
< HitPairGeneratorFromLayerPair
thePairGenerator
 

Detailed Description

A HitTripletGenerator from HitPairGenerator and vector of Layers. The HitPairGenerator provides a set of hit pairs. For each pair the search for compatible hit(s) is done among provided Layers

Definition at line 19 of file PixelTripletLargeTipGenerator.h.

Member Typedef Documentation

Definition at line 21 of file PixelTripletLargeTipGenerator.h.

Constructor & Destructor Documentation

PixelTripletLargeTipGenerator::PixelTripletLargeTipGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)

Definition at line 48 of file PixelTripletLargeTipGenerator.cc.

50  useFixedPreFiltering(cfg.getParameter<bool>("useFixedPreFiltering")),
51  extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
52  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
53  useMScat(cfg.getParameter<bool>("useMultScattering")),
54  useBend(cfg.getParameter<bool>("useBending")),
55  dphi(useFixedPreFiltering ? cfg.getParameter<double>("phiPreFiltering") : 0) {}
T getParameter(std::string const &) const
PixelTripletLargeTipGenerator::~PixelTripletLargeTipGenerator ( )
virtual

Definition at line 57 of file PixelTripletLargeTipGenerator.cc.

57 {}

Member Function Documentation

void PixelTripletLargeTipGenerator::hitTriplets ( const TrackingRegion region,
OrderedHitTriplets trs,
const edm::Event ev,
const edm::EventSetup es,
SeedingLayerSetsHits::SeedingLayerSet  pairLayers,
const std::vector< SeedingLayerSetsHits::SeedingLayer > &  thirdLayers 
)
overridevirtual

Implements HitTripletGeneratorFromPairAndLayers.

Definition at line 73 of file PixelTripletLargeTipGenerator.cc.

References angle(), checkPhiInRange(), PixelRecoUtilities::curvature(), declareDynArray, dphi, relativeConstraints::empty, ev, extraHitRPhitolerance, extraHitRZtolerance, f, plotBeamSpotDB::first, fnSigmaRZ, Geom::ftwoPi(), edm::EventSetup::get(), i, HitDoublets::inner, DetLayer::isBarrel(), geometryCSVtoXML::line, M_PI, bookConverter::max, ThirdHitPredictionFromCircle::HelixRZ::maxCurvature(), min(), normalizedPhi(), nSigmaPhi, nSigmaRZ, PixelRecoLineRZ::origin(), TrackingRegion::originRBound(), HitDoublets::outer, p3, PV3DBase< T, PVType, FrameType >::perp(), GeometricSearchDet::position(), edm::ESHandle< class >::product(), proxim(), TrackingRegion::ptMin(), Basic2DVector< T >::r(), CosmicsPD_Skims::radius, HLT_25ns10e33_v2_cff::region, OrderedHitTriplets::size(), findQualityFiles::size, std::swap(), HitTripletGeneratorFromPairAndLayers::theMaxElement, HitTripletGeneratorFromPairAndLayers::thePairGenerator, patCandidatesForDimuonsSequences_cff::tracker, useBend, useFixedPreFiltering, useMScat, findQualityFiles::v, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

79 {
81  es.get<TrackerDigiGeometryRecord>().get(tracker);
82 
83  //Retrieve tracker topology from geometry
85  es.get<TrackerTopologyRcd>().get(tTopoHand);
86  const TrackerTopology *tTopo=tTopoHand.product();
87 
88  auto const & doublets = thePairGenerator->doublets(region,ev,es, pairLayers);
89 
90  if (doublets.empty()) return;
91 
92  auto outSeq = doublets.detLayer(HitDoublets::outer)->seqNum();
93 
94 
95  int size = thirdLayers.size();
96 
97 
98  using NodeInfo = KDTreeNodeInfo<unsigned int>;
99  std::vector<NodeInfo > layerTree; // re-used throughout
100  std::vector<unsigned int> foundNodes; // re-used throughout
101  foundNodes.reserve(100);
102 
104  declareDynArray(LayerRZPredictions, size, mapPred);
105 
106  float rzError[size]; //save maximum errors
107 
108  const float maxDelphi = region.ptMin() < 0.3f ? float(M_PI)/4.f : float(M_PI)/8.f; // FIXME move to config??
109  const float maxphi = M_PI+maxDelphi, minphi = -maxphi; // increase to cater for any range
110  const float safePhi = M_PI-maxDelphi; // sideband
111 
112 
113  const RecHitsSortedInPhi * thirdHitMap[size];
114 
115  for(int il = 0; il < size; il++) {
116  thirdHitMap[il] = &(*theLayerCache)(thirdLayers[il], region, ev, es);
117  auto const & hits = *thirdHitMap[il];
118 
119  const DetLayer *layer = thirdLayers[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  predRZ.correction.init(es, region.ptMin(), *doublets.detLayer(HitDoublets::inner), *doublets.detLayer(HitDoublets::outer), *thirdLayers[il].detLayer(), useMScat, false);
129 
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  // populate side-bands
143  if (angle>safePhi) layerTree.emplace_back(i, angle-Geom::ftwoPi(), v);
144  else if (angle<-safePhi) layerTree.emplace_back(i, angle+Geom::ftwoPi(), v);
145  }
146  KDTreeBox phiZ(minphi, maxphi, minv-0.01f, maxv+0.01f); // declare our bounds
147  //add fudge factors in case only one hit and also for floating-point inaccuracy
148  hitTree[il].build(layerTree, phiZ); // make KDtree
149  rzError[il] = maxErr; //save error
150  }
151 
152  double curv = PixelRecoUtilities::curvature(1. / region.ptMin(), es);
153 
154  for (std::size_t ip =0; ip!=doublets.size(); ip++) {
155  auto xi = doublets.x(ip,HitDoublets::inner);
156  auto yi = doublets.y(ip,HitDoublets::inner);
157  auto zi = doublets.z(ip,HitDoublets::inner);
158  // auto rvi = doublets.rv(ip,HitDoublets::inner);
159  auto xo = doublets.x(ip,HitDoublets::outer);
160  auto yo = doublets.y(ip,HitDoublets::outer);
161  auto zo = doublets.z(ip,HitDoublets::outer);
162  // auto rvo = doublets.rv(ip,HitDoublets::outer);
163  GlobalPoint gp1(xi,yi,zi);
164  GlobalPoint gp2(xo,yo,zo);
165 
166  auto toPos = std::signbit(zo-zi);
167 
168  PixelRecoLineRZ line(gp1, gp2);
169  PixelRecoPointRZ point2(gp2.perp(), zo);
170  ThirdHitPredictionFromCircle predictionRPhi(gp1, gp2, extraHitRPhitolerance);
171 
172  Range generalCurvature = predictionRPhi.curvature(region.originRBound());
173  if (!intersect(generalCurvature, Range(-curv, curv))) continue;
174 
175  for(int il = 0; il < size; il++) {
176  if (hitTree[il].empty()) continue; // Don't bother if no hits
177  const DetLayer *layer = thirdLayers[il].detLayer();
178  bool barrelLayer = layer->isBarrel();
179 
180  if ( (!barrelLayer) & (toPos != std::signbit(layer->position().z())) ) continue;
181 
182 
183  Range curvature = generalCurvature;
184 
185  LayerRZPredictions &predRZ = mapPred[il];
186  predRZ.line.initPropagator(&line);
187 
188  auto & correction = predRZ.correction;
189  correction.init(line, point2, outSeq);
190 
191 
192  Range rzRange;
193  if (useBend) {
194  // For the barrel region:
195  // swiping the helix passing through the two points across from
196  // negative to positive bending, can give us a sort of U-shaped
197  // projection onto the phi-z (barrel) or r-z plane (forward)
198  // so we checking minimum/maximum of all three possible extrema
199  //
200  // For the endcap region:
201  // Checking minimum/maximum radius of the helix projection
202  // onto an endcap plane, here we have to guard against
203  // looping tracks, when phi(delta z) gets out of control.
204  // HelixRZ::rAtZ should not follow looping tracks, but clamp
205  // to the minimum reachable r with the next-best lower |curvature|.
206  // So same procedure as for the barrel region can be applied.
207  //
208  // In order to avoid looking for potential looping tracks at all
209  // we also clamp the allowed curvature range for this layer,
210  // and potentially fail the layer entirely
211 
212  if (!barrelLayer) {
213  Range z3s = predRZ.line.detRange();
214  double z3 = z3s.first < 0 ? std::max(z3s.first, z3s.second)
215  : std::min(z3s.first, z3s.second);
216  double maxCurvature = HelixRZ::maxCurvature(&predictionRPhi,
217  gp1.z(), gp2.z(), z3);
218  if (!intersect(curvature, Range(-maxCurvature, maxCurvature)))
219  continue;
220  }
221 
222  HelixRZ helix1(&predictionRPhi, gp1.z(), gp2.z(), curvature.first);
223  HelixRZ helix2(&predictionRPhi, gp1.z(), gp2.z(), curvature.second);
224 
225  predRZ.helix1.initPropagator(&helix1);
226  predRZ.helix2.initPropagator(&helix2);
227 
228  Range rzRanges[2] = { predRZ.helix1(), predRZ.helix2() };
229  predRZ.helix1.initPropagator(nullptr);
230  predRZ.helix2.initPropagator(nullptr);
231 
232  rzRange.first = std::min(rzRanges[0].first, rzRanges[1].first);
233  rzRange.second = std::max(rzRanges[0].second, rzRanges[1].second);
234 
235  // if the allowed curvatures include a straight line,
236  // this can give us another extremum for allowed r/z
237  if (curvature.first * curvature.second < 0.0) {
238  Range rzLineRange = predRZ.line();
239  rzRange.first = std::min(rzRange.first, rzLineRange.first);
240  rzRange.second = std::max(rzRange.second, rzLineRange.second);
241  }
242  } else {
243  rzRange = predRZ.line();
244  }
245 
246  if (rzRange.first >= rzRange.second)
247  continue;
248 
249  correction.correctRZRange(rzRange);
250 
251  Range phiRange;
252  if (useFixedPreFiltering) {
253  float phi0 = doublets.phi(ip,HitDoublets::outer);
254  phiRange = Range(phi0 - dphi, phi0 + dphi);
255  } else {
256  Range radius;
257 
258  if (barrelLayer) {
259  radius = predRZ.line.detRange();
260  if (!intersect(rzRange, predRZ.line.detSize()))
261  continue;
262  } else {
263  radius = rzRange;
264  if (!intersect(radius, predRZ.line.detSize()))
265  continue;
266  }
267 
268  //gc: predictionRPhi uses the cosine rule to find the phi of the 3rd point at radius, assuming the pairCurvature range [-c,+c]
269  if ( (curvature.first<0.0f) & (curvature.second<0.0f) ) {
270  radius.swap();
271  } else if ( (curvature.first>=0.0f) & (curvature.second>=0.0f) ) {;}
272  else {
273  radius.first=radius.second;
274  }
275  auto phi12 = predictionRPhi.phi(curvature.first,radius.first);
276  auto phi22 = predictionRPhi.phi(curvature.second,radius.second);
277  phi12 = normalizedPhi(phi12);
278  phi22 = proxim(phi22,phi12);
279  phiRange = Range(phi12,phi22); phiRange.sort();
280  auto rmean = radius.mean();
281  phiRange.first *= rmean;
282  phiRange.second *= rmean;
283  correction.correctRPhiRange(phiRange);
284  phiRange.first /= rmean;
285  phiRange.second /= rmean;
286 
287  }
288 
289  foundNodes.clear(); // Now recover hits in bounding box...
290  float prmin=phiRange.min(), prmax=phiRange.max(); //get contiguous range
291 
292  if (prmax-prmin>maxDelphi) {
293  auto prm = phiRange.mean();
294  prmin = prm - 0.5f*maxDelphi;
295  prmax = prm + 0.5f*maxDelphi;
296  }
297 
298  if (barrelLayer) {
299  Range regMax = predRZ.line.detRange();
300  Range regMin = predRZ.line(regMax.min());
301  regMax = predRZ.line(regMax.max());
302  correction.correctRZRange(regMin);
303  correction.correctRZRange(regMax);
304  if (regMax.min() < regMin.min()) { std::swap(regMax, regMin);}
305  KDTreeBox phiZ(prmin, prmax,
306  regMin.min()-fnSigmaRZ*rzError[il],
307  regMax.max()+fnSigmaRZ*rzError[il]);
308  hitTree[il].search(phiZ, foundNodes);
309  }
310  else {
311  KDTreeBox phiZ(prmin, prmax,
312  rzRange.min()-fnSigmaRZ*rzError[il],
313  rzRange.max()+fnSigmaRZ*rzError[il]);
314  hitTree[il].search(phiZ, foundNodes);
315  }
316 
317  MatchedHitRZCorrectionFromBending l2rzFixup(doublets.hit(ip,HitDoublets::outer)->det()->geographicalId(), tTopo);
318  MatchedHitRZCorrectionFromBending l3rzFixup = predRZ.rzPositionFixup;
319 
320  thirdHitMap[il] = &(*theLayerCache)(thirdLayers[il], region, ev, es);
321  auto const & hits = *thirdHitMap[il];
322  for (auto KDdata : foundNodes) {
323  GlobalPoint p3 = hits.gp(KDdata);
324  double p3_r = p3.perp();
325  double p3_z = p3.z();
326  float p3_phi = hits.phi(KDdata);
327 
328  Range rangeRPhi = predictionRPhi(curvature, p3_r);
329  correction.correctRPhiRange(rangeRPhi);
330 
331  float ir = 1.f/p3_r;
332  float phiErr = nSigmaPhi * hits.drphi[KDdata]*ir;
333  if (!checkPhiInRange(p3_phi, rangeRPhi.first*ir-phiErr, rangeRPhi.second*ir+phiErr))
334  continue;
335 
336  Basic2DVector<double> thc(p3.x(), p3.y());
337 
338  auto curv_ = predictionRPhi.curvature(thc);
339  double p2_r = point2.r(); double p2_z = point2.z(); // they will be modified!
340 
341  l2rzFixup(predictionRPhi, curv_, *doublets.hit(ip,HitDoublets::outer), p2_r, p2_z, tTopo);
342  l3rzFixup(predictionRPhi, curv_, *hits.theHits[KDdata].hit(), p3_r, p3_z, tTopo);
343 
344  Range rangeRZ;
345  if (useBend) {
346  HelixRZ updatedHelix(&predictionRPhi, gp1.z(), p2_z, curv_);
347  rangeRZ = predRZ.helix1(barrelLayer ? p3_r : p3_z, updatedHelix);
348  } else {
349  float tIP = predictionRPhi.transverseIP(thc);
350  PixelRecoPointRZ updatedPoint2(p2_r, p2_z);
351  PixelRecoLineRZ updatedLine(line.origin(), point2, tIP);
352  rangeRZ = predRZ.line(barrelLayer ? p3_r : p3_z, line);
353  }
354  correction.correctRZRange(rangeRZ);
355 
356  double err = nSigmaRZ * hits.dv[KDdata];
357 
358  rangeRZ.first -= err, rangeRZ.second += err;
359 
360  if (!rangeRZ.inside(barrelLayer ? p3_z : p3_r)) continue;
361 
362  if (theMaxElement!=0 && result.size() >= theMaxElement) {
363  result.clear();
364  edm::LogError("TooManyTriplets")<<" number of triples exceed maximum. no triplets produced.";
365  return;
366  }
367  result.emplace_back( doublets.hit(ip,HitDoublets::inner), doublets.hit(ip,HitDoublets::outer), hits.theHits[KDdata].hit());
368  }
369  }
370  }
371  // std::cout << "found triplets " << result.size() << std::endl;
372 }
int i
Definition: DBlmapReader.cc:9
T perp() const
Definition: PV3DBase.h:72
constexpr float ftwoPi()
Definition: Pi.h:36
PixelRecoRange< float > Range
T y() const
Definition: PV3DBase.h:63
T r() const
Radius, same as mag()
bool ev
tuple result
Definition: mps_fire.py:83
U second(std::pair< T, U > const &p)
T curvature(T InversePt, const edm::EventSetup &iSetup)
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
T z() const
Definition: PV3DBase.h:64
double f[11][100]
T min(T a, T b)
Definition: MathUtil.h:58
#define M_PI
bool isBarrel() const
Definition: DetLayer.h:32
static double maxCurvature(const ThirdHitPredictionFromCircle *circle, double z1, double z2, double z3)
std::unique_ptr< HitPairGeneratorFromLayerPair > thePairGenerator
T proxim(T b, T a)
Definition: normalizedPhi.h:13
virtual const Surface::PositionType & position() const
Returns position of the surface.
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
float ptMin() const
minimal pt of interest
bool checkPhiInRange(T phi, T phi1, T phi2)
Definition: normalizedPhi.h:21
#define declareDynArray(T, n, x)
Definition: DynArray.h:58
T x() const
Definition: PV3DBase.h:62
tuple size
Write out results.
T normalizedPhi(T phi)
Definition: normalizedPhi.h:8
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

Member Data Documentation

const float PixelTripletLargeTipGenerator::dphi
private

Definition at line 39 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

const float PixelTripletLargeTipGenerator::extraHitRPhitolerance
private

Definition at line 36 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

const float PixelTripletLargeTipGenerator::extraHitRZtolerance
private

Definition at line 35 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

const bool PixelTripletLargeTipGenerator::useBend
private

Definition at line 38 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

const bool PixelTripletLargeTipGenerator::useFixedPreFiltering
private

Definition at line 34 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().

const bool PixelTripletLargeTipGenerator::useMScat
private

Definition at line 37 of file PixelTripletLargeTipGenerator.h.

Referenced by hitTriplets().