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 Member Functions | Private Attributes
PixelTripletHLTGenerator Class Reference

#include <PixelTripletHLTGenerator.h>

Inheritance diagram for PixelTripletHLTGenerator:
HitTripletGeneratorFromPairAndLayers HitTripletGenerator OrderedHitsGenerator

Public Member Functions

virtual void hitTriplets (const TrackingRegion &region, OrderedHitTriplets &trs, const edm::Event &ev, const edm::EventSetup &es)
 
virtual void init (const HitPairGenerator &pairs, const std::vector< ctfseeding::SeedingLayer > &layers, LayerCacheType *layerCache)
 
const HitPairGeneratorpairGenerator () const
 
 PixelTripletHLTGenerator (const edm::ParameterSet &cfg)
 
const std::vector
< ctfseeding::SeedingLayer > & 
thirdLayers () const
 
virtual ~PixelTripletHLTGenerator ()
 
- Public Member Functions inherited from HitTripletGeneratorFromPairAndLayers
virtual ~HitTripletGeneratorFromPairAndLayers ()
 
- Public Member Functions inherited from HitTripletGenerator
virtual void clear ()
 
 HitTripletGenerator (unsigned int size=500)
 
virtual void hitTriplets (const TrackingRegion &reg, OrderedHitTriplets &prs, const edm::EventSetup &es)
 
virtual const OrderedHitTripletsrun (const TrackingRegion &region, const edm::Event &ev, const edm::EventSetup &es)
 
virtual ~HitTripletGenerator ()
 
- Public Member Functions inherited from OrderedHitsGenerator
 OrderedHitsGenerator ()
 
virtual ~OrderedHitsGenerator ()
 

Private Types

typedef
CombinedHitTripletGenerator::LayerCacheType 
LayerCacheType
 

Private Member Functions

bool checkPhiInRange (float phi, float phi1, float phi2) const
 
std::pair< float, float > mergePhiRanges (const std::pair< float, float > &r1, const std::pair< float, float > &r2) const
 

Private Attributes

float dphi
 
float extraHitRPhitolerance
 
float extraHitRZtolerance
 
SeedComparitortheComparitor
 
LayerCacheTypetheLayerCache
 
std::vector
< ctfseeding::SeedingLayer
theLayers
 
HitPairGeneratorthePairGenerator
 
bool useBend
 
bool useFixedPreFiltering
 
bool useMScat
 

Additional Inherited Members

- Public Types inherited from HitTripletGeneratorFromPairAndLayers
typedef LayerHitMapCache LayerCacheType
 
- Public Attributes inherited from OrderedHitsGenerator
unsigned int theMaxElement
 

Detailed Description

Definition at line 23 of file PixelTripletHLTGenerator.h.

Member Typedef Documentation

Definition at line 25 of file PixelTripletHLTGenerator.h.

Constructor & Destructor Documentation

PixelTripletHLTGenerator::PixelTripletHLTGenerator ( const edm::ParameterSet cfg)

Definition at line 32 of file PixelTripletHLTGenerator.cc.

References dphi, reco::get(), edm::ParameterSet::getParameter(), AlCaHLTBitMon_QueryRunRegistry::string, theComparitor, OrderedHitsGenerator::theMaxElement, and useFixedPreFiltering.

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 }
T getParameter(std::string const &) const
T get(const Candidate &c)
Definition: component.h:55
PixelTripletHLTGenerator::~PixelTripletHLTGenerator ( )
virtual

Definition at line 51 of file PixelTripletHLTGenerator.cc.

References theComparitor, and thePairGenerator.

51  {
52  delete thePairGenerator;
53  delete theComparitor;
54 }

Member Function Documentation

bool PixelTripletHLTGenerator::checkPhiInRange ( float  phi,
float  phi1,
float  phi2 
) const
private

Definition at line 305 of file PixelTripletHLTGenerator.cc.

References Geom::ftwoPi().

Referenced by hitTriplets().

306 {
307  while (phi > phi2) phi -= Geom::ftwoPi();
308  while (phi < phi1) phi += Geom::ftwoPi();
309  return ( (phi1 <= phi) && (phi <= phi2) );
310 }
float ftwoPi()
Definition: Pi.h:36
Definition: DDAxes.h:10
void PixelTripletHLTGenerator::hitTriplets ( const TrackingRegion region,
OrderedHitTriplets trs,
const edm::Event ev,
const edm::EventSetup es 
)
virtual

Implements HitTripletGenerator.

Definition at line 65 of file PixelTripletHLTGenerator.cc.

References angle(), KDTreeLinkerAlgo< DATA >::build(), checkPhiInRange(), SeedComparitor::compatible(), constexpr, PixelRecoUtilities::curvature(), HitDoublets::detLayer(), HitPairGenerator::doublets(), dphi, relativeConstraints::empty, extraHitRPhitolerance, extraHitRZtolerance, f, Geom::fpi(), Geom::ftwoPi(), i, SeedComparitor::init(), ThirdHitRZPredictionBase::initLayer(), ThirdHitRZPrediction< Propagator >::initPropagator(), ThirdHitRZPredictionBase::initTolerance(), HitDoublets::inner, PixelRecoRange< T >::intersection(), DetLayer::isBarrel(), geometryCSVtoXML::line, LogDebug, max(), mergePhiRanges(), bookConverter::min, nSigmaPhi, nSigmaRZ, TrackingRegion::origin(), TrackingRegion::originRBound(), HitDoublets::outer, PV3DBase< T, PVType, FrameType >::perp(), TrackingRegion::ptMin(), CosmicsPD_Skims::radius, KDTreeLinkerAlgo< DATA >::search(), DetLayer::seqNum(), OrderedHitTriplets::size(), findQualityFiles::size, mathSSE::sqrt(), swap(), theComparitor, theLayers, OrderedHitsGenerator::theMaxElement, thePairGenerator, useBend, useFixedPreFiltering, useMScat, findQualityFiles::v, PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

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  std::vector<unsigned int> foundNodes; // re-used thoughout
93  foundNodes.reserve(100);
94 
96  float rzError[size]; //save maximum errors
97  float maxphi = Geom::ftwoPi(), minphi = -maxphi; // increase to cater for any range
98 
99  // fill the prediction vector
100  for (int il=0; il!=size; ++il) {
101  thirdHitMap[il] = &(*theLayerCache)(&theLayers[il], region, ev, es);
102  auto const & hits = *thirdHitMap[il];
103  ThirdHitRZPrediction<PixelRecoLineRZ> & pred = preds[il];
104  pred.initLayer(theLayers[il].detLayer());
106 
107  layerTree.clear();
108  float minv=999999.0, maxv= -999999.0; // Initialise to extreme values in case no hits
109  float maxErr=0.0f;
110  for (unsigned int i=0; i!=hits.size(); ++i) {
111  auto angle = hits.phi(i);
112  auto v = hits.gv(i);
113  //use (phi,r) for endcaps rather than (phi,z)
114  minv = std::min(minv,v); maxv = std::max(maxv,v);
115  float myerr = hits.dv[i];
116  maxErr = std::max(maxErr,myerr);
117  layerTree.emplace_back(i, angle, v); // save it
118  if (angle < 0) // wrap all points in phi
119  { layerTree.emplace_back(i, angle+Geom::ftwoPi(), v);}
120  else
121  { layerTree.emplace_back(i, angle-Geom::ftwoPi(), v);}
122  }
123  KDTreeBox phiZ(minphi, maxphi, minv-0.01f, maxv+0.01f); // declare our bounds
124  //add fudge factors in case only one hit and also for floating-point inaccuracy
125  hitTree[il].build(layerTree, phiZ); // make KDtree
126  rzError[il] = maxErr; //save error
127  // std::cout << "layer " << theLayers[il].detLayer()->seqNum() << " " << layerTree.size() << std::endl;
128  }
129 
130  float imppar = region.originRBound();
131  float imppartmp = region.originRBound()+region.origin().perp();
132  float curv = PixelRecoUtilities::curvature(1.f/region.ptMin(), es);
133 
134  for (std::size_t ip =0; ip!=doublets.size(); ip++) {
135  auto xi = doublets.x(ip,HitDoublets::inner);
136  auto yi = doublets.y(ip,HitDoublets::inner);
137  auto zi = doublets.z(ip,HitDoublets::inner);
138  auto rvi = doublets.rv(ip,HitDoublets::inner);
139  auto xo = doublets.x(ip,HitDoublets::outer);
140  auto yo = doublets.y(ip,HitDoublets::outer);
141  auto zo = doublets.z(ip,HitDoublets::outer);
142  auto rvo = doublets.rv(ip,HitDoublets::outer);
143 
144  PixelRecoPointRZ point1(rvi, zi);
145  PixelRecoPointRZ point2(rvo, zo);
146  PixelRecoLineRZ line(point1, point2);
147  ThirdHitPredictionFromInvParabola predictionRPhi(xi-region.origin().x(),yi-region.origin().y(),
148  xo-region.origin().x(),yo-region.origin().y(),
149  imppar,curv,extraHitRPhitolerance);
150 
151  ThirdHitPredictionFromInvParabola predictionRPhitmp(xi,yi,xo,yo,imppartmp,curv,extraHitRPhitolerance);
152 
153  // printf("++Constr %f %f %f %f %f %f %f\n",xi,yi,xo,yo,imppartmp,curv,extraHitRPhitolerance);
154 
155  // std::cout << ip << ": " << point1.r() << ","<< point1.z() << " "
156  // << point2.r() << ","<< point2.z() <<std::endl;
157 
158  for (int il=0; il!=size; ++il) {
159  if (hitTree[il].empty()) continue; // Don't bother if no hits
160 
161  auto const & hits = *thirdHitMap[il];
162 
163  const DetLayer * layer = theLayers[il].detLayer();
164  auto barrelLayer = layer->isBarrel();
165 
166  ThirdHitCorrection correction(es, region.ptMin(), layer, line, point2, outSeq, useMScat, useBend);
167 
168  ThirdHitRZPrediction<PixelRecoLineRZ> & predictionRZ = preds[il];
169 
170  predictionRZ.initPropagator(&line);
171  Range rzRange = predictionRZ();
172  correction.correctRZRange(rzRange);
173 
174  Range phiRange;
175  if (useFixedPreFiltering) {
176  float phi0 = doublets.phi(ip,HitDoublets::outer);
177  phiRange = Range(phi0-dphi,phi0+dphi);
178  }
179  else {
180  Range radius;
181  if (barrelLayer) {
182  radius = predictionRZ.detRange();
183  } else {
184  radius = Range(max(rzRange.min(), predictionRZ.detSize().min()),
185  min(rzRange.max(), predictionRZ.detSize().max()) );
186  }
187  if (radius.empty()) continue;
188 
189  // std::cout << "++R " << radius.min() << " " << radius.max() << std::endl;
190 
191  /*
192  Range rPhi1m = predictionRPhitmp(radius.max(), -1);
193  Range rPhi1p = predictionRPhitmp(radius.max(), 1);
194  Range rPhi2m = predictionRPhitmp(radius.min(), -1);
195  Range rPhi2p = predictionRPhitmp(radius.min(), 1);
196  Range rPhi1 = rPhi1m.sum(rPhi1p);
197  Range rPhi2 = rPhi2m.sum(rPhi2p);
198 
199 
200  auto rPhi1N = predictionRPhitmp(radius.max());
201  auto rPhi2N = predictionRPhitmp(radius.min());
202 
203  std::cout << "VI "
204  << rPhi1N.first <<'/'<< rPhi1.first << ' '
205  << rPhi1N.second <<'/'<< rPhi1.second << ' '
206  << rPhi2N.first <<'/'<< rPhi2.first << ' '
207  << rPhi2N.second <<'/'<< rPhi2.second
208  << std::endl;
209 
210  */
211 
212  auto rPhi1 = predictionRPhitmp(radius.max());
213  auto rPhi2 = predictionRPhitmp(radius.min());
214 
215 
216  correction.correctRPhiRange(rPhi1);
217  correction.correctRPhiRange(rPhi2);
218  rPhi1.first /= radius.max();
219  rPhi1.second /= radius.max();
220  rPhi2.first /= radius.min();
221  rPhi2.second /= radius.min();
222  phiRange = mergePhiRanges(rPhi1,rPhi2);
223  }
224 
225  constexpr float nSigmaRZ = std::sqrt(12.f); // ...and continue as before
226  constexpr float nSigmaPhi = 3.f;
227 
228  foundNodes.clear(); // Now recover hits in bounding box...
229  float prmin=phiRange.min(), prmax=phiRange.max();
230  if ((prmax-prmin) > Geom::ftwoPi())
231  { prmax=Geom::fpi(); prmin = -Geom::fpi();}
232  else
233  { while (prmax>maxphi) { prmin -= Geom::ftwoPi(); prmax -= Geom::ftwoPi();}
234  while (prmin<minphi) { prmin += Geom::ftwoPi(); prmax += Geom::ftwoPi();}
235  // This needs range -twoPi to +twoPi to work
236  }
237  if (barrelLayer)
238  {
239  Range regMax = predictionRZ.detRange();
240  Range regMin = predictionRZ(regMax.min()-regOffset);
241  regMax = predictionRZ(regMax.max()+regOffset);
242  correction.correctRZRange(regMin);
243  correction.correctRZRange(regMax);
244  if (regMax.min() < regMin.min()) { swap(regMax, regMin);}
245  KDTreeBox phiZ(prmin, prmax, regMin.min()-nSigmaRZ*rzError[il], regMax.max()+nSigmaRZ*rzError[il]);
246  hitTree[il].search(phiZ, foundNodes);
247  }
248  else
249  {
250  KDTreeBox phiZ(prmin, prmax,
251  rzRange.min()-regOffset-nSigmaRZ*rzError[il],
252  rzRange.max()+regOffset+nSigmaRZ*rzError[il]);
253  hitTree[il].search(phiZ, foundNodes);
254  }
255 
256  // std::cout << ip << ": " << theLayers[il].detLayer()->seqNum() << " " << foundNodes.size() << " " << prmin << " " << prmax << std::endl;
257 
258 
259  // int kk=0;
260  for (auto KDdata : foundNodes) {
261 
262  if (theMaxElement!=0 && result.size() >= theMaxElement){
263  result.clear();
264  edm::LogError("TooManyTriplets")<<" number of triples exceeds maximum. no triplets produced.";
265  return;
266  }
267 
268  float p3_u = hits.u[KDdata];
269  float p3_v = hits.v[KDdata];
270  float p3_phi = hits.lphi[KDdata];
271 
272  //if ((kk++)%100==0)
273  //std::cout << kk << ": " << p3_u << " " << p3_v << " " << p3_phi << std::endl;
274 
275 
276  Range allowed = predictionRZ(p3_u);
277  correction.correctRZRange(allowed);
278  float vErr = nSigmaRZ *hits.dv[KDdata];
279  Range hitRange(p3_v-vErr, p3_v+vErr);
280  Range crossingRange = allowed.intersection(hitRange);
281  if (crossingRange.empty()) continue;
282 
283  float ir = 1.f/hits.rv(KDdata);
284  float phiErr = nSigmaPhi * hits.drphi[KDdata]*ir;
285  for (int icharge=-1; icharge <=1; icharge+=2) {
286  Range rangeRPhi = predictionRPhi(hits.rv(KDdata), icharge);
287  correction.correctRPhiRange(rangeRPhi);
288  if (checkPhiInRange(p3_phi, rangeRPhi.first*ir-phiErr, rangeRPhi.second*ir+phiErr)) {
289  // insert here check with comparitor
290  OrderedHitTriplet hittriplet( doublets.hit(ip,HitDoublets::inner), doublets.hit(ip,HitDoublets::outer), hits.theHits[KDdata].hit());
291  if (!theComparitor || theComparitor->compatible(hittriplet,region) ) {
292  result.push_back( hittriplet );
293  } else {
294  LogDebug("RejectedTriplet") << "rejected triplet from comparitor ";
295  }
296  break;
297  }
298  }
299  }
300  }
301  }
302  // std::cout << "triplets " << result.size() << std::endl;
303 }
#define LogDebug(id)
void swap(ora::Record &rh, ora::Record &lh)
Definition: Record.h:70
int i
Definition: DBlmapReader.cc:9
std::vector< ctfseeding::SeedingLayer > theLayers
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
void initLayer(const DetLayer *layer)
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
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
PixelRecoRange< T > intersection(const PixelRecoRange< T > &r) const
bool checkPhiInRange(float phi, float phi1, float phi2) const
TransientTrackingRecHit::ConstRecHitPointer Hit
void initTolerance(float tolerance)
DetLayer const * detLayer(layer l) const
tuple size
Write out results.
#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
void PixelTripletHLTGenerator::init ( const HitPairGenerator pairs,
const std::vector< ctfseeding::SeedingLayer > &  layers,
LayerCacheType layerCache 
)
virtual

Implements HitTripletGeneratorFromPairAndLayers.

Definition at line 56 of file PixelTripletHLTGenerator.cc.

References HitPairGenerator::clone(), theLayerCache, theLayers, and thePairGenerator.

59 {
60  thePairGenerator = pairs.clone();
61  theLayers = layers;
62  theLayerCache = layerCache;
63 }
std::vector< ctfseeding::SeedingLayer > theLayers
virtual HitPairGenerator * clone() const =0
std::pair< float, float > PixelTripletHLTGenerator::mergePhiRanges ( const std::pair< float, float > &  r1,
const std::pair< float, float > &  r2 
) const
private

Definition at line 312 of file PixelTripletHLTGenerator.cc.

References Geom::fpi(), Geom::ftwoPi(), max(), and bookConverter::min.

Referenced by hitTriplets().

314 { float r2_min=r2.first;
315  float r2_max=r2.second;
316  while (r1.first-r2_min > Geom::fpi()) { r2_min += Geom::ftwoPi(); r2_max += Geom::ftwoPi();}
317  while (r1.first-r2_min < -Geom::fpi()) { r2_min -= Geom::ftwoPi(); r2_max -= Geom::ftwoPi(); }
318 
319  return std::make_pair(min(r1.first,r2_min),max(r1.second,r2_max));
320 }
float fpi()
Definition: Pi.h:35
const T & max(const T &a, const T &b)
float ftwoPi()
Definition: Pi.h:36
const HitPairGenerator& PixelTripletHLTGenerator::pairGenerator ( ) const
inline

Definition at line 38 of file PixelTripletHLTGenerator.h.

References thePairGenerator.

38 { return *thePairGenerator; }
const std::vector<ctfseeding::SeedingLayer>& PixelTripletHLTGenerator::thirdLayers ( ) const
inline

Definition at line 39 of file PixelTripletHLTGenerator.h.

References theLayers.

39 { return theLayers; }
std::vector< ctfseeding::SeedingLayer > theLayers

Member Data Documentation

float PixelTripletHLTGenerator::dphi
private

Definition at line 56 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets(), and PixelTripletHLTGenerator().

float PixelTripletHLTGenerator::extraHitRPhitolerance
private

Definition at line 53 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

float PixelTripletHLTGenerator::extraHitRZtolerance
private

Definition at line 52 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

SeedComparitor* PixelTripletHLTGenerator::theComparitor
private
LayerCacheType* PixelTripletHLTGenerator::theLayerCache
private

Definition at line 49 of file PixelTripletHLTGenerator.h.

Referenced by init().

std::vector<ctfseeding::SeedingLayer> PixelTripletHLTGenerator::theLayers
private

Definition at line 48 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets(), init(), and thirdLayers().

HitPairGenerator* PixelTripletHLTGenerator::thePairGenerator
private
bool PixelTripletHLTGenerator::useBend
private

Definition at line 55 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

bool PixelTripletHLTGenerator::useFixedPreFiltering
private

Definition at line 51 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets(), and PixelTripletHLTGenerator().

bool PixelTripletHLTGenerator::useMScat
private

Definition at line 54 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().