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