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 }
Vector3DBase
Definition: Vector3DBase.h:8
HLT_2018_cff.extraHitRPhitolerance
extraHitRPhitolerance
Definition: HLT_2018_cff.py:8543
PixelTripletNoTipGenerator::extraHitRZtolerance
float extraHitRZtolerance
Definition: PixelTripletNoTipGenerator.h:38
OrderedHitPairs
Definition: OrderedHitPairs.h:8
ZMuMuCategoriesSequences_cff.chi2Cut
chi2Cut
Definition: ZMuMuCategoriesSequences_cff.py:141
photonAnalyzer_cfi.zMax
zMax
Definition: photonAnalyzer_cfi.py:95
DetLayer
Definition: DetLayer.h:21
BarrelDetLayer::surface
const BoundSurface & surface() const final
GeometricSearchDet interface.
Definition: BarrelDetLayer.h:29
ForwardDetLayer::surface
const BoundSurface & surface() const final
The surface of the GeometricSearchDet.
Definition: ForwardDetLayer.h:29
edm::swap
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:117
PixelRecoUtilities::inversePt
T inversePt(T curvature, const edm::EventSetup &iSetup)
Definition: PixelRecoUtilities.h:48
PixelTripletNoTipGenerator_cfi.extraHitPhiToleranceForPreFiltering
extraHitPhiToleranceForPreFiltering
Definition: PixelTripletNoTipGenerator_cfi.py:6
GeometricSearchDet::position
virtual const Surface::PositionType & position() const
Returns position of the surface.
Definition: GeometricSearchDet.h:31
PixelRecoUtilities::curvature
T curvature(T InversePt, const edm::EventSetup &iSetup)
Definition: PixelRecoUtilities.h:42
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
GeomDetEnumerators::barrel
Definition: GeomDetEnumerators.h:9
DetLayer::location
virtual Location location() const =0
Which part of the detector (barrel, endcap)
HLT_2018_cff.doublets
doublets
Definition: HLT_2018_cff.py:8544
RecHitsSortedInPhi.h
createJobs.theNSigma
theNSigma
Definition: createJobs.py:313
sqr
T sqr(T t)
Definition: PixelTripletNoTipGenerator.cc:21
beampixel_dqm_sourceclient-live_cfg.zRange
zRange
Definition: beampixel_dqm_sourceclient-live_cfg.py:148
HitTripletGeneratorFromPairAndLayers::thePairGenerator
std::unique_ptr< HitPairGeneratorFromLayerPair > thePairGenerator
Definition: HitTripletGeneratorFromPairAndLayers.h:55
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
Surface::bounds
const Bounds & bounds() const
Definition: Surface.h:87
HitTripletGeneratorFromPairAndLayers
Definition: HitTripletGeneratorFromPairAndLayers.h:25
p2
double p2[4]
Definition: TauolaWrapper.h:90
PixelTripletNoTipGenerator::hitTriplets
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
Definition: PixelTripletNoTipGenerator.cc:37
RecHitsSortedInPhi::Hit
BaseTrackerRecHit const * Hit
Definition: RecHitsSortedInPhi.h:19
GlobalPoint
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
PixelTripletNoTipGenerator::PixelTripletNoTipGenerator
PixelTripletNoTipGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
Definition: PixelTripletNoTipGenerator.cc:27
PixelRecoUtilities.h
Point3DBase< float, GlobalTag >
OrderedSet.t
t
Definition: OrderedSet.py:90
RecHitsSortedInPhi
Definition: RecHitsSortedInPhi.h:17
PixelRecoRange< float >
Bounds::thickness
virtual float thickness() const =0
PixelTripletNoTipGenerator::theChi2Cut
double theChi2Cut
Definition: PixelTripletNoTipGenerator.h:42
HitDoublets
Definition: RecHitsSortedInPhi.h:124
edm::ParameterSet
Definition: ParameterSet.h:36
edm::LogError
Definition: MessageLogger.h:183
Event.h
jetfilter_cfi.ptCut
ptCut
Definition: jetfilter_cfi.py:6
OrderedHitTriplets
Definition: OrderedHitTriplets.h:9
ThirdHitZPrediction.h
PixelTripletNoTipGenerator::~PixelTripletNoTipGenerator
~PixelTripletNoTipGenerator() override
Definition: PixelTripletNoTipGenerator.cc:35
HLT_2018_cff.extraHitRZtolerance
extraHitRZtolerance
Definition: HLT_2018_cff.py:44971
photonAnalyzer_cfi.zMin
zMin
Definition: photonAnalyzer_cfi.py:94
PixelRecoLineRZ
Definition: PixelRecoLineRZ.h:12
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:50
BarrelDetLayer.h
p1
double p1[4]
Definition: TauolaWrapper.h:89
edm::EventSetup
Definition: EventSetup.h:57
BarrelDetLayer
Definition: BarrelDetLayer.h:22
ThirdHitZPrediction
Definition: ThirdHitZPrediction.h:8
PixelRecoPointRZ
Definition: PixelRecoPointRZ.h:6
looper.cfg
cfg
Definition: looper.py:297
MultipleScatteringParametrisation.h
ThirdHitPredictionFromInvLine.h
std
Definition: JetResolutionObject.h:76
PixelTripletNoTipGenerator.h
RecHitsSortedInPhi::hits
std::vector< Hit > hits(float phiMin, float phiMax) const
Definition: RecHitsSortedInPhi.cc:93
HitTripletGeneratorFromPairAndLayers::theMaxElement
const unsigned int theMaxElement
Definition: HitTripletGeneratorFromPairAndLayers.h:57
ForwardDetLayer
Definition: ForwardDetLayer.h:22
HitPairGeneratorFromLayerPair.h
edm::shift
static unsigned const int shift
Definition: LuminosityBlockID.cc:7
T
long double T
Definition: Basic3DVectorLD.h:48
ev
bool ev
Definition: Hydjet2Hadronizer.cc:95
ForwardDetLayer.h
flavorHistoryFilter_cfi.dr
dr
Definition: flavorHistoryFilter_cfi.py:37
Exception
Definition: hltDiff.cc:246
SeedingLayerSetsHits::SeedingLayerSet
Definition: SeedingLayerSetsHits.h:65
ThirdHitPredictionFromInvLine
Definition: ThirdHitPredictionFromInvLine.h:18
CosmicsPD_Skims.radius
radius
Definition: CosmicsPD_Skims.py:135
PixelRecoRange::sort
PixelRecoRange< T > & sort()
Definition: PixelRecoRange.h:51
HLT_2018_cff.region
region
Definition: HLT_2018_cff.py:81479
Range
PixelRecoRange< float > Range
Definition: PixelTripletNoTipGenerator.cc:19
DetLayer.h
OrderedHitTriplet
Definition: OrderedHitTriplet.h:11
EventSetup.h
p3
double p3[4]
Definition: TauolaWrapper.h:91
TrackingRegion
Definition: TrackingRegion.h:40
MultipleScatteringParametrisation
Definition: MultipleScatteringParametrisation.h:16
mps_fire.result
result
Definition: mps_fire.py:303
SimpleDiskBounds.h
ParameterSet.h
edm::Event
Definition: Event.h:73
mps_splice.line
line
Definition: mps_splice.py:76
BarrelDetLayer::specificSurface
virtual const BoundCylinder & specificSurface() const final
Extension of the interface.
Definition: BarrelDetLayer.h:39
PixelTripletNoTipGenerator::extraHitPhiToleranceForPreFiltering
float extraHitPhiToleranceForPreFiltering
Definition: PixelTripletNoTipGenerator.h:40
PV3DBase::perp
T perp() const
Definition: PV3DBase.h:69
edm::ConsumesCollector
Definition: ConsumesCollector.h:39
PV3DBase::phi
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
Hit
SeedingHitSet::ConstRecHitPointer Hit
Definition: SeedGeneratorFromProtoTracksEDProducer.cc:34
findQualityFiles.size
size
Write out results.
Definition: findQualityFiles.py:443
PixelTripletNoTipGenerator::theNSigma
double theNSigma
Definition: PixelTripletNoTipGenerator.h:41