CMS 3D CMS Logo

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