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 
7 
10 #include "ThirdHitCorrection.h"
13 #include <iostream>
14 
17 
19 
22 
23 using namespace std;
24 using namespace ctfseeding;
25 
27  : thePairGenerator(0),
28  theLayerCache(0),
29  useFixedPreFiltering(cfg.getParameter<bool>("useFixedPreFiltering")),
30  extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
31  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
32  useMScat(cfg.getParameter<bool>("useMultScattering")),
33  useBend(cfg.getParameter<bool>("useBending"))
34 {
35  theMaxElement=cfg.getParameter<unsigned int>("maxElement");
36  dphi = (useFixedPreFiltering) ? cfg.getParameter<double>("phiPreFiltering") : 0;
37 
38  edm::ParameterSet comparitorPSet =
39  cfg.getParameter<edm::ParameterSet>("SeedComparitorPSet");
40  std::string comparitorName = comparitorPSet.getParameter<std::string>("ComponentName");
41  theComparitor = (comparitorName == "none") ?
42  0 : SeedComparitorFactory::get()->create( comparitorName, comparitorPSet);
43 
44 }
45 
47 
48 { delete thePairGenerator;
49  delete theComparitor;
50 }
51 
53  const std::vector<SeedingLayer> & layers,
54  LayerCacheType* layerCache)
55 {
56  thePairGenerator = pairs.clone();
57  theLayers = layers;
58  theLayerCache = layerCache;
59 }
60 
62  const TrackingRegion& region,
64  const edm::Event & ev,
65  const edm::EventSetup& es)
66 {
68  OrderedHitPairs pairs; pairs.reserve(30000);
69  OrderedHitPairs::const_iterator ip;
70 
71  thePairGenerator->hitPairs(region,pairs,ev,es);
72 
73  if (pairs.empty()) return;
74 
75  int size = theLayers.size();
76 
77  typedef std::vector<ThirdHitRZPrediction<PixelRecoLineRZ> > Preds;
78  Preds preds(size);
79 
80  std::vector<const RecHitsSortedInPhi *> thirdHitMap(size);
82  vector<Hit> thirdHits;
83 
84  // fill the prediciton vetor
85  for (int il=0; il!=size; ++il) {
86  thirdHitMap[il] = &(*theLayerCache)(&theLayers[il], region, ev, es);
87  ThirdHitRZPrediction<PixelRecoLineRZ> & pred = preds[il];
88  pred.initLayer(theLayers[il].detLayer());
89  pred.initTolerance(extraHitRZtolerance);
90  }
91 
92 
93  double imppar = region.originRBound();
94  double curv = PixelRecoUtilities::curvature(1/region.ptMin(), es);
95 
96  for (ip = pairs.begin(); ip != pairs.end(); ip++) {
97 
98  GlobalPoint gp1tmp = (*ip).inner()->globalPosition();
99  GlobalPoint gp2tmp = (*ip).outer()->globalPosition();
100  GlobalPoint gp1(gp1tmp.x()-region.origin().x(), gp1tmp.y()-region.origin().y(), gp1tmp.z());
101  GlobalPoint gp2(gp2tmp.x()-region.origin().x(), gp2tmp.y()-region.origin().y(), gp2tmp.z());
102 
103  PixelRecoPointRZ point1(gp1.perp(), gp1.z());
104  PixelRecoPointRZ point2(gp2.perp(), gp2.z());
105  PixelRecoLineRZ line(point1, point2);
106  ThirdHitPredictionFromInvParabola predictionRPhi(gp1,gp2,imppar,curv,extraHitRPhitolerance);
107  ThirdHitPredictionFromInvParabola predictionRPhitmp(gp1tmp,gp2tmp,imppar+region.origin().perp(),curv,extraHitRPhitolerance);
108 
109 
110  for (int il=0; il!=size; ++il) {
111  const DetLayer * layer = theLayers[il].detLayer();
112 // bool pixelLayer = ( layer->subDetector() == GeomDetEnumerators::PixelBarrel
113 // || layer->subDetector() == GeomDetEnumerators::PixelEndcap);
114  bool barrelLayer = (layer->location() == GeomDetEnumerators::barrel);
115 
116  ThirdHitCorrection correction(es, region.ptMin(), layer, line, point2, useMScat, useBend);
117 
118  ThirdHitRZPrediction<PixelRecoLineRZ> & predictionRZ = preds[il];
119 
120  predictionRZ.initPropagator(&line);
121  Range rzRange = predictionRZ();
122 
123  correction.correctRZRange(rzRange);
124  Range phiRange;
125  if (useFixedPreFiltering) {
126  float phi0 = (*ip).outer()->globalPosition().phi();
127  phiRange = Range(phi0-dphi,phi0+dphi);
128  }
129  else {
130  Range radius;
131  if (barrelLayer) {
132  radius = predictionRZ.detRange();
133  } else {
134  radius = Range(
135  max(rzRange.min(), predictionRZ.detSize().min()),
136  min(rzRange.max(), predictionRZ.detSize().max()) );
137  }
138  if (radius.empty()) continue;
139  Range rPhi1m = predictionRPhitmp(radius.max(), -1);
140  Range rPhi1p = predictionRPhitmp(radius.max(), 1);
141  Range rPhi2m = predictionRPhitmp(radius.min(), -1);
142  Range rPhi2p = predictionRPhitmp(radius.min(), 1);
143  Range rPhi1 = rPhi1m.sum(rPhi1p);
144  Range rPhi2 = rPhi2m.sum(rPhi2p);
145  correction.correctRPhiRange(rPhi1);
146  correction.correctRPhiRange(rPhi2);
147  rPhi1.first /= radius.max();
148  rPhi1.second /= radius.max();
149  rPhi2.first /= radius.min();
150  rPhi2.second /= radius.min();
151  phiRange = mergePhiRanges(rPhi1,rPhi2);
152  }
153 
154 // LayerHitMapLoop thirdHits =
155 // pixelLayer ? thirdHitMap[il]->loop(phiRange, rzRange) :
156 // thirdHitMap[il]->loop();
157 
158  thirdHits.clear();
159  thirdHitMap[il]->hits(phiRange.min(),phiRange.max(), thirdHits);
160 
161  static float nSigmaRZ = std::sqrt(12.f);
162  static float nSigmaPhi = 3.f;
163 
164  typedef vector<Hit>::const_iterator IH;
165  for (IH th=thirdHits.begin(), eh=thirdHits.end(); th !=eh; ++th) {
166 
167  if (theMaxElement!=0 && result.size() >= theMaxElement){
168  result.clear();
169  edm::LogError("TooManyTriplets")<<" number of triples exceed maximum. no triplets produced.";
170  return;
171  }
172  const Hit& hit = (*th);
173  GlobalPoint point(hit->globalPosition().x()-region.origin().x(),
174  hit->globalPosition().y()-region.origin().y(),
175  hit->globalPosition().z() );
176  float p3_r = point.perp();
177  float p3_z = point.z();
178  float p3_phi = point.phi();
179 
180  if (barrelLayer) {
181  Range allowedZ = predictionRZ(p3_r);
182  correction.correctRZRange(allowedZ);
183 
184  float zErr = nSigmaRZ * hit->errorGlobalZ();
185  Range hitRange(p3_z-zErr, p3_z+zErr);
186  Range crossingRange = allowedZ.intersection(hitRange);
187  if (crossingRange.empty()) continue;
188  } else {
189  Range allowedR = predictionRZ(p3_z);
190  correction.correctRZRange(allowedR);
191  float rErr = nSigmaRZ * hit->errorGlobalR();
192  Range hitRange(p3_r-rErr, p3_r+rErr);
193  Range crossingRange = allowedR.intersection(hitRange);
194  if (crossingRange.empty()) continue;
195  }
196 
197  float phiErr = nSigmaPhi*hit->errorGlobalRPhi()/p3_r;
198  for (int icharge=-1; icharge <=1; icharge+=2) {
199  Range rangeRPhi = predictionRPhi(p3_r, icharge);
200  correction.correctRPhiRange(rangeRPhi);
201  if (checkPhiInRange(p3_phi, rangeRPhi.first/p3_r-phiErr, rangeRPhi.second/p3_r+phiErr)) {
202  // insert here check with comparitor
203  OrderedHitTriplet hittriplet( (*ip).inner(), (*ip).outer(), hit);
204  if(!theComparitor || theComparitor->compatible(hittriplet,region) ) {
205  result.push_back( hittriplet );
206  } else {
207  LogDebug("RejectedTriplet") << "rejected triplet from comparitor "
208  << hittriplet.outer()->globalPosition().x() << " "
209  << hittriplet.outer()->globalPosition().y() << " "
210  << hittriplet.outer()->globalPosition().z();
211  }
212  break;
213  }
214  }
215  }
216  }
217  }
218 
219 }
220 
221 bool PixelTripletHLTGenerator::checkPhiInRange(float phi, float phi1, float phi2) const
222 {
223  while (phi > phi2) phi -= Geom::ftwoPi();
224  while (phi < phi1) phi += Geom::ftwoPi();
225  return ( (phi1 <= phi) && (phi <= phi2) );
226 }
227 
229  const std::pair<float,float>& r1, const std::pair<float,float>& r2) const
230 {
231  float r2_min=r2.first;
232  float r2_max=r2.second;
233  while (r1.first-r2_min > Geom::fpi()) { r2_min += Geom::ftwoPi(); r2_max += Geom::ftwoPi();}
234  while (r1.first-r2_min < -Geom::fpi()) { r2_min -= Geom::ftwoPi(); r2_max -= Geom::ftwoPi(); }
235 
236  return std::make_pair(min(r1.first,r2_min),max(r1.second,r2_max));
237 }
#define LogDebug(id)
T getParameter(std::string const &) const
std::vector< ctfseeding::SeedingLayer > theLayers
virtual HitPairGenerator * clone() const =0
void initPropagator(const Propagator *propagator)
T perp() const
Definition: PV3DBase.h:71
virtual float ptMin() const =0
minimal pt of interest
virtual GlobalPoint origin() const =0
T max() const
void initLayer(const DetLayer *layer)
virtual Location location() const =0
Which part of the detector (barrel, endcap)
static const double nSigmaPhi
PixelRecoRange< T > sum(const PixelRecoRange< T > &r) const
T y() const
Definition: PV3DBase.h:62
#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)
std::pair< float, float > mergePhiRanges(const std::pair< float, float > &r1, const std::pair< float, float > &r2) const
T curvature(T InversePt, const edm::EventSetup &iSetup)
const T & max(const T &a, const T &b)
virtual void hitPairs(const TrackingRegion &reg, OrderedHitPairs &prs, const edm::EventSetup &es)
T sqrt(T t)
Definition: SSEVec.h:46
T z() const
Definition: PV3DBase.h:63
tuple result
Definition: query.py:137
double f[11][100]
static const double nSigmaRZ
PixelRecoRange< float > Range
PixelTripletHLTGenerator(const edm::ParameterSet &cfg)
virtual float originRBound() const =0
bounds the particle vertex in the transverse plane
PixelRecoRange< T > intersection(const PixelRecoRange< T > &r) const
bool checkPhiInRange(float phi, float phi1, float phi2) const
TransientTrackingRecHit::ConstRecHitPointer Hit
virtual unsigned int size() const
T x() const
Definition: PV3DBase.h:61
virtual void hitTriplets(const TrackingRegion &region, OrderedHitTriplets &trs, const edm::Event &ev, const edm::EventSetup &es)
tuple size
Write out results.
T get(const Candidate &c)
Definition: component.h:56
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
float ftwoPi()
Definition: Pi.h:36
Definition: DDAxes.h:10