CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelTripletNoTipGenerator.cc
Go to the documentation of this file.
4 
9 #include "ThirdHitZPrediction.h"
10 
13 
18 
20 template<class T> T sqr(T t) { return t * t;}
21 
22 using namespace std;
23 using namespace ctfseeding;
24 
26  : thePairGenerator(0),
27  theLayerCache(0),
28  extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
29  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
30  extraHitPhiToleranceForPreFiltering(cfg.getParameter<double>("extraHitPhiToleranceForPreFiltering")),
31  theNSigma(cfg.getParameter<double>("nSigma")),
32  theChi2Cut(cfg.getParameter<double>("chi2Cut"))
33 { }
34 
36  const std::vector<SeedingLayer> & layers,
37  LayerCacheType* layerCache)
38 {
39  thePairGenerator = pairs.clone();
40  theLayers = layers;
41  theLayerCache = layerCache;
42 }
43 
45  const TrackingRegion& region,
47  const edm::Event & ev,
48  const edm::EventSetup& es)
49 {
50 
51 //
52 //edm::Handle<reco::BeamSpot> bsHandle;
53 //ev.getByLabel( theBeamSpotTag, bsHandle);
54 //if(!bsHandle.isValid()) return;
55 //const reco::BeamSpot & bs = *bsHandle;
56 //double errorXY = sqrt( sqr(bs.BeamWidthX()) + sqr(bs.BeamWidthY()) );
57 //
58 
59  GlobalPoint bsPoint = region.origin();
60  double errorXY = region.originRBound();
61  GlobalVector shift = bsPoint - GlobalPoint(0.,0.,0.);
62 
63  OrderedHitPairs pairs; pairs.reserve(30000);
64  OrderedHitPairs::const_iterator ip;
65  thePairGenerator->hitPairs(region,pairs,ev,es);
66 
67  if (pairs.size() ==0) return;
68 
69  int size = theLayers.size();
70 
71  const RecHitsSortedInPhi **thirdHitMap = new const RecHitsSortedInPhi*[size];
72  for (int il=0; il <=size-1; il++) {
73  thirdHitMap[il] = &(*theLayerCache)(&theLayers[il], region, ev, es);
74  }
75 
76  const HitPairGeneratorFromLayerPair * pairGen = dynamic_cast<const HitPairGeneratorFromLayerPair *>(thePairGenerator);
77  const DetLayer * firstLayer = pairGen->innerLayer().detLayer();
78  const DetLayer * secondLayer = pairGen->outerLayer().detLayer();
79  if (!firstLayer || !secondLayer) return;
80 
81  MultipleScatteringParametrisation sigma1RPhi( firstLayer, es);
82  MultipleScatteringParametrisation sigma2RPhi( secondLayer, es);
83 
85  for (ip = pairs.begin(); ip != pairs.end(); ip++) {
86 
87  GlobalPoint p1((*ip).inner()->globalPosition()-shift);
88  GlobalPoint p2((*ip).outer()->globalPosition()-shift);
89 
90  ThirdHitPredictionFromInvLine predictionRPhiTMP(p1, p2 );
91  double pt_p1p2 = 1./PixelRecoUtilities::inversePt(predictionRPhiTMP.curvature(),es);
92 
93  PixelRecoPointRZ point1(p1.perp(), p1.z());
94  PixelRecoPointRZ point2(p2.perp(), p2.z());
95 
96  PixelRecoLineRZ line(point1, point2);
97  double msRPhi1 = sigma1RPhi(pt_p1p2, line.cotLine());
98  double msRPhi2 = sigma2RPhi(pt_p1p2, line.cotLine(),point1);
99  double sinTheta = 1/sqrt(1+sqr(line.cotLine()));
100  double cosTheta = fabs(line.cotLine())/sqrt(1+sqr(line.cotLine()));
101 
102  double p1_errorRPhi = sqrt(sqr((*ip).inner()->errorGlobalRPhi())+sqr(msRPhi1) +sqr(errorXY));
103  double p2_errorRPhi = sqrt(sqr((*ip).outer()->errorGlobalRPhi())+sqr(msRPhi2) +sqr(errorXY));
104 
105  ThirdHitPredictionFromInvLine predictionRPhi(p1, p2, p1_errorRPhi, p2_errorRPhi );
106 
107  for (int il=0; il <=size-1; il++) {
108 
109  const DetLayer * layer = theLayers[il].detLayer();
110  bool barrelLayer = (layer->location() == GeomDetEnumerators::barrel);
111  MultipleScatteringParametrisation sigma3RPhi( layer, es);
112  double msRPhi3 = sigma3RPhi(pt_p1p2, line.cotLine(),point2);
113 
114  Range rRange;
115  if (barrelLayer) {
116  const BarrelDetLayer& bl = dynamic_cast<const BarrelDetLayer&>(*layer);
117  float halfThickness = bl.surface().bounds().thickness()/2;
118  float radius = bl.specificSurface().radius();
119  rRange = Range(radius-halfThickness, radius+halfThickness);
120  } else {
121  const ForwardDetLayer& fl = dynamic_cast<const ForwardDetLayer&>(*layer);
122  float halfThickness = fl.surface().bounds().thickness()/2;
123  float zLayer = fl.position().z() ;
124  float zMin = zLayer-halfThickness;
125  float zMax = zLayer+halfThickness;
126  GlobalVector dr = p2-p1;
127  GlobalPoint p3_a = p2+dr*(zMin-p2.z())/dr.z();
128  GlobalPoint p3_b = p2+dr*(zMax-p2.z())/dr.z();
129  if (zLayer * p3_a.z() < 0) continue;
130  rRange = Range(p3_a.perp(), p3_b.perp());
131  rRange.sort();
132  }
133  double displacment = shift.perp();
134  GlobalPoint crossing1 = predictionRPhi.crossing(rRange.min()-displacment)+shift;
135  GlobalPoint crossing2 = predictionRPhi.crossing(rRange.max()+displacment)+shift;
136  float c1_phi= crossing1.phi();
137  float c2_phi= crossing2.phi();
138  if (c2_phi < c1_phi) swap(c1_phi,c2_phi);
139  if (c2_phi-c1_phi > M_PI) { c2_phi -= 2*M_PI; swap(c1_phi,c2_phi); }
140  double extraAngle = (displacment+theNSigma*msRPhi3)/rRange.min()+extraHitPhiToleranceForPreFiltering;
141  c1_phi -= extraAngle;
142  c2_phi += extraAngle;
143  vector<Hit> thirdHits = thirdHitMap[il]->hits(c1_phi, c2_phi) ;
144 
145  typedef vector<Hit>::const_iterator IH;
146  for (IH th=thirdHits.begin(), eh=thirdHits.end(); th < eh; ++th) {
147  GlobalPoint p3((*th)->globalPosition()-shift);
148  double p3_errorRPhi = sqrt(sqr((*th)->errorGlobalRPhi()) +sqr(msRPhi3) + sqr(errorXY));
149 
150  predictionRPhi.add(p3,p3_errorRPhi);
151 
152  double curvature = predictionRPhi.curvature();
153 
154  ThirdHitZPrediction zPrediction( (*ip).inner()->globalPosition(), sqrt(sqr((*ip).inner()->errorGlobalR())+sqr(msRPhi1/cosTheta)), sqrt( sqr((*ip).inner()->errorGlobalZ())+ sqr(msRPhi1/sinTheta)),
155  (*ip).outer()->globalPosition(), sqrt(sqr((*ip).outer()->errorGlobalR())+sqr(msRPhi2/cosTheta)), sqrt( sqr((*ip).outer()->errorGlobalZ())+sqr(msRPhi2/sinTheta)),
156  curvature, theNSigma);
157  ThirdHitZPrediction::Range zRange = zPrediction((*th)->globalPosition(), sqrt(sqr((*th)->errorGlobalR()))+sqr(msRPhi3/cosTheta));
158 
159  double z3Hit = (*th)->globalPosition().z();
160  double z3HitError = theNSigma*(sqrt(sqr((*th)->errorGlobalZ()) + sqr(msRPhi3/sinTheta) ))+extraHitRZtolerance;
161  Range hitZRange(z3Hit-z3HitError, z3Hit+z3HitError);
162  bool inside = hitZRange.hasIntersection(zRange);
163 
164  double curvatureMS = PixelRecoUtilities::curvature(1./region.ptMin(),es);
165  bool ptCut = (predictionRPhi.curvature()-theNSigma*predictionRPhi.errorCurvature() < curvatureMS);
166  bool chi2Cut = (predictionRPhi.chi2() < theChi2Cut);
167  if (inside && ptCut && chi2Cut) {
168  if (theMaxElement!=0 && result.size() >= theMaxElement){
169  result.clear();
170  edm::LogError("TooManyTriplets")<<" number of triples exceed maximum. no triplets produced.";
171  delete [] thirdHitMap;
172  return;
173  }
174  result.push_back( OrderedHitTriplet( (*ip).inner(), (*ip).outer(), *th));
175  }
176  predictionRPhi.remove(p3,p3_errorRPhi);
177  }
178  }
179  }
180  delete [] thirdHitMap;
181 }
void swap(ora::Record &rh, ora::Record &lh)
Definition: Record.h:70
virtual void hitTriplets(const TrackingRegion &region, OrderedHitTriplets &trs, const edm::Event &ev, const edm::EventSetup &es)
virtual HitPairGenerator * clone() const =0
PixelTripletNoTipGenerator(const edm::ParameterSet &cfg)
T perp() const
Definition: PV3DBase.h:71
virtual float ptMin() const =0
minimal pt of interest
virtual GlobalPoint origin() const =0
virtual Location location() const =0
Which part of the detector (barrel, endcap)
const DetLayer * detLayer() const
Definition: SeedingLayer.cc:80
std::vector< Hit > hits(float phiMin, float phiMax) const
virtual const BoundSurface & surface() const
The surface of the GeometricSearchDet.
Geom::Phi< T > phi() const
Definition: PV3DBase.h:68
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::vector< ctfseeding::SeedingLayer > theLayers
virtual unsigned int size() const
virtual float thickness() const =0
T inversePt(T curvature, const edm::EventSetup &iSetup)
T curvature(T InversePt, const edm::EventSetup &iSetup)
Scalar radius() const
Radius of the cylinder.
Definition: Cylinder.h:55
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
virtual void init(const HitPairGenerator &pairs, const std::vector< ctfseeding::SeedingLayer > &layers, LayerCacheType *layerCache)
double p2[4]
Definition: TauolaWrapper.h:90
const Bounds & bounds() const
Definition: BoundSurface.h:89
PixelRecoRange< float > Range
#define M_PI
Definition: BFit3D.cc:3
virtual const BoundSurface & surface() const
GeometricSearchDet interface.
virtual const Surface::PositionType & position() const
Returns position of the surface.
virtual const BoundCylinder & specificSurface() const
Extension of the interface.
virtual float originRBound() const =0
bounds the particle vertex in the transverse plane
double p1[4]
Definition: TauolaWrapper.h:89
Square< F >::type sqr(const F &f)
Definition: Square.h:13
TransientTrackingRecHit::ConstRecHitPointer Hit
static unsigned int const shift
long double T
virtual unsigned int size() const
tuple size
Write out results.
double p3[4]
Definition: TauolaWrapper.h:91