CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Attributes
PixelTripletNoTipGenerator Class Reference

#include <PixelTripletNoTipGenerator.h>

Inheritance diagram for PixelTripletNoTipGenerator:
HitTripletGeneratorFromPairAndLayers

Public Member Functions

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
 
void hitTriplets (const TrackingRegion &region, OrderedHitTriplets &result, const edm::EventSetup &es, const HitDoublets &doublets, const RecHitsSortedInPhi **thirdHitMap, const std::vector< const DetLayer * > &thirdLayerDetLayer, const int nThirdLayers) override
 
 PixelTripletNoTipGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
 
 ~PixelTripletNoTipGenerator () override
 
- Public Member Functions inherited from HitTripletGeneratorFromPairAndLayers
 HitTripletGeneratorFromPairAndLayers (unsigned int maxElement=0)
 
 HitTripletGeneratorFromPairAndLayers (const edm::ParameterSet &pset)
 
void init (std::unique_ptr< HitPairGeneratorFromLayerPair > &&pairs, LayerCacheType *layerCache)
 
const HitPairGeneratorFromLayerPairpairGenerator () const
 
virtual ~HitTripletGeneratorFromPairAndLayers ()
 

Private Types

typedef CombinedHitTripletGenerator::LayerCacheType LayerCacheType
 

Private Attributes

float extraHitPhiToleranceForPreFiltering
 
float extraHitRPhitolerance
 
float extraHitRZtolerance
 
double theChi2Cut
 
double theNSigma
 

Additional Inherited Members

- Public Types inherited from HitTripletGeneratorFromPairAndLayers
typedef LayerHitMapCache LayerCacheType
 
- Static Public Member Functions inherited from HitTripletGeneratorFromPairAndLayers
static void fillDescriptions (edm::ParameterSetDescription &desc)
 
- Protected Attributes inherited from HitTripletGeneratorFromPairAndLayers
LayerCacheTypetheLayerCache
 
const unsigned int theMaxElement
 
std::unique_ptr< HitPairGeneratorFromLayerPairthePairGenerator
 

Detailed Description

Definition at line 13 of file PixelTripletNoTipGenerator.h.

Member Typedef Documentation

Definition at line 14 of file PixelTripletNoTipGenerator.h.

Constructor & Destructor Documentation

PixelTripletNoTipGenerator::PixelTripletNoTipGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)

Definition at line 24 of file PixelTripletNoTipGenerator.cc.

26  extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
27  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
28  extraHitPhiToleranceForPreFiltering(cfg.getParameter<double>("extraHitPhiToleranceForPreFiltering")),
29  theNSigma(cfg.getParameter<double>("nSigma")),
30  theChi2Cut(cfg.getParameter<double>("chi2Cut"))
31 { }
T getParameter(std::string const &) const
PixelTripletNoTipGenerator::~PixelTripletNoTipGenerator ( )
override

Definition at line 33 of file PixelTripletNoTipGenerator.cc.

33 {}

Member Function Documentation

void PixelTripletNoTipGenerator::hitTriplets ( const TrackingRegion region,
OrderedHitTriplets trs,
const edm::Event ev,
const edm::EventSetup es,
const SeedingLayerSetsHits::SeedingLayerSet pairLayers,
const std::vector< SeedingLayerSetsHits::SeedingLayer > &  thirdLayers 
)
overridevirtual

Implements HitTripletGeneratorFromPairAndLayers.

Definition at line 35 of file PixelTripletNoTipGenerator.cc.

References GeomDetEnumerators::barrel, Surface::bounds(), PixelTripletNoTipGenerator_cfi::chi2Cut, PixelRecoUtilities::curvature(), runTauDisplay::dr, extraHitPhiToleranceForPreFiltering, extraHitRZtolerance, RecHitsSortedInPhi::hits(), PixelRecoUtilities::inversePt(), mps_splice::line, DetLayer::location(), M_PI, TrackingRegion::origin(), TrackingRegion::originRBound(), p1, p2, p3, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), GeometricSearchDet::position(), HiGenCleaner_cff::ptCut, TrackingRegion::ptMin(), TCMET_cfi::radius, edm::shift, OrderedHitTriplets::size(), findQualityFiles::size, PixelRecoRange< T >::sort(), BarrelDetLayer::specificSurface(), sqr(), mathSSE::sqrt(), ForwardDetLayer::surface(), BarrelDetLayer::surface(), edm::swap(), theChi2Cut, HitTripletGeneratorFromPairAndLayers::theMaxElement, theNSigma, HitTripletGeneratorFromPairAndLayers::thePairGenerator, Bounds::thickness(), PV3DBase< T, PVType, FrameType >::z(), hfnoseDigiStudy_cfi::zMax, and hfnoseDigiStudy_cfi::zMin.

42 {
43 
44 //
45 //edm::Handle<reco::BeamSpot> bsHandle;
46 //ev.getByLabel( theBeamSpotTag, bsHandle);
47 //if(!bsHandle.isValid()) return;
48 //const reco::BeamSpot & bs = *bsHandle;
49 //double errorXY = sqrt( sqr(bs.BeamWidthX()) + sqr(bs.BeamWidthY()) );
50 //
51 
52  const GlobalPoint& bsPoint = region.origin();
53  double errorXY = region.originRBound();
54  GlobalVector shift = bsPoint - GlobalPoint(0.,0.,0.);
55 
56  OrderedHitPairs pairs; pairs.reserve(30000);
57  OrderedHitPairs::const_iterator ip;
58  thePairGenerator->hitPairs(region,pairs,ev,es, pairLayers);
59 
60  if (pairs.empty()) return;
61 
62  int size = thirdLayers.size();
63 
64  const RecHitsSortedInPhi **thirdHitMap = new const RecHitsSortedInPhi*[size];
65  for (int il=0; il <=size-1; il++) {
66  thirdHitMap[il] = &(*theLayerCache)(thirdLayers[il], region, es);
67  }
68 
69  const DetLayer * firstLayer = thePairGenerator->innerLayer(pairLayers).detLayer();
70  const DetLayer * secondLayer = thePairGenerator->outerLayer(pairLayers).detLayer();
71  if (!firstLayer || !secondLayer) return;
72 
73  MultipleScatteringParametrisation sigma1RPhi( firstLayer, es);
74  MultipleScatteringParametrisation sigma2RPhi( secondLayer, es);
75 
77  for (ip = pairs.begin(); ip != pairs.end(); ip++) {
78 
79  GlobalPoint p1((*ip).inner()->globalPosition()-shift);
80  GlobalPoint p2((*ip).outer()->globalPosition()-shift);
81 
82  ThirdHitPredictionFromInvLine predictionRPhiTMP(p1, p2 );
83  double pt_p1p2 = 1./PixelRecoUtilities::inversePt(predictionRPhiTMP.curvature(),es);
84 
85  PixelRecoPointRZ point1(p1.perp(), p1.z());
86  PixelRecoPointRZ point2(p2.perp(), p2.z());
87 
88  PixelRecoLineRZ line(point1, point2);
89  double msRPhi1 = sigma1RPhi(pt_p1p2, line.cotLine(), 0.f);
90  double msRPhi2 = sigma2RPhi(pt_p1p2, line.cotLine(),point1);
91  double sinTheta = 1/sqrt(1+sqr(line.cotLine()));
92  double cosTheta = fabs(line.cotLine())/sqrt(1+sqr(line.cotLine()));
93 
94  double p1_errorRPhi = sqrt(sqr((*ip).inner()->errorGlobalRPhi())+sqr(msRPhi1) +sqr(errorXY));
95  double p2_errorRPhi = sqrt(sqr((*ip).outer()->errorGlobalRPhi())+sqr(msRPhi2) +sqr(errorXY));
96 
97  ThirdHitPredictionFromInvLine predictionRPhi(p1, p2, p1_errorRPhi, p2_errorRPhi );
98 
99  for (int il=0; il <=size-1; il++) {
100 
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) continue;
122  rRange = Range(p3_a.perp(), p3_b.perp());
123  rRange.sort();
124  }
125  double displacment = shift.perp();
126  GlobalPoint crossing1 = predictionRPhi.crossing(rRange.min()-displacment)+shift;
127  GlobalPoint crossing2 = predictionRPhi.crossing(rRange.max()+displacment)+shift;
128  float c1_phi= crossing1.phi();
129  float c2_phi= crossing2.phi();
130  if (c2_phi < c1_phi) swap(c1_phi,c2_phi);
131  if (c2_phi-c1_phi > M_PI) { c2_phi -= 2*M_PI; swap(c1_phi,c2_phi); }
132  double extraAngle = (displacment+theNSigma*msRPhi3)/rRange.min()+extraHitPhiToleranceForPreFiltering;
133  c1_phi -= extraAngle;
134  c2_phi += extraAngle;
135  vector<Hit> thirdHits = thirdHitMap[il]->hits(c1_phi, c2_phi) ;
136 
137  typedef vector<Hit>::const_iterator IH;
138  for (IH th=thirdHits.begin(), eh=thirdHits.end(); th < eh; ++th) {
139  GlobalPoint p3((*th)->globalPosition()-shift);
140  double p3_errorRPhi = sqrt(sqr((*th)->errorGlobalRPhi()) +sqr(msRPhi3) + sqr(errorXY));
141 
142  predictionRPhi.add(p3,p3_errorRPhi);
143 
144  double curvature = predictionRPhi.curvature();
145 
146  ThirdHitZPrediction zPrediction( (*ip).inner()->globalPosition(), sqrt(sqr((*ip).inner()->errorGlobalR())+sqr(msRPhi1/cosTheta)), sqrt( sqr((*ip).inner()->errorGlobalZ())+ sqr(msRPhi1/sinTheta)),
147  (*ip).outer()->globalPosition(), sqrt(sqr((*ip).outer()->errorGlobalR())+sqr(msRPhi2/cosTheta)), sqrt( sqr((*ip).outer()->errorGlobalZ())+sqr(msRPhi2/sinTheta)),
148  curvature, theNSigma);
149  ThirdHitZPrediction::Range zRange = zPrediction((*th)->globalPosition(), sqrt(sqr((*th)->errorGlobalR()))+sqr(msRPhi3/cosTheta));
150 
151  double z3Hit = (*th)->globalPosition().z();
152  double z3HitError = theNSigma*(sqrt(sqr((*th)->errorGlobalZ()) + sqr(msRPhi3/sinTheta) ))+extraHitRZtolerance;
153  Range hitZRange(z3Hit-z3HitError, z3Hit+z3HitError);
154  bool inside = hitZRange.hasIntersection(zRange);
155 
156  double curvatureMS = PixelRecoUtilities::curvature(1./region.ptMin(),es);
157  bool ptCut = (predictionRPhi.curvature()-theNSigma*predictionRPhi.errorCurvature() < curvatureMS);
158  bool chi2Cut = (predictionRPhi.chi2() < theChi2Cut);
159  if (inside && ptCut && chi2Cut) {
160  if (theMaxElement!=0 && result.size() >= theMaxElement){
161  result.clear();
162  edm::LogError("TooManyTriplets")<<" number of triples exceed maximum. no triplets produced.";
163  delete [] thirdHitMap;
164  return;
165  }
166  result.push_back( OrderedHitTriplet( (*ip).inner(), (*ip).outer(), *th));
167  }
168  predictionRPhi.remove(p3,p3_errorRPhi);
169  }
170  }
171  }
172  delete [] thirdHitMap;
173 }
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:72
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:69
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
const Bounds & bounds() const
Definition: Surface.h:120
virtual const BoundCylinder & specificSurface() const final
Extension of the interface.
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:116
T inversePt(T curvature, const edm::EventSetup &iSetup)
T curvature(T InversePt, const edm::EventSetup &iSetup)
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
double p2[4]
Definition: TauolaWrapper.h:90
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
const BoundSurface & surface() const final
GeometricSearchDet interface.
double p3[4]
Definition: TauolaWrapper.h:91
void PixelTripletNoTipGenerator::hitTriplets ( const TrackingRegion region,
OrderedHitTriplets result,
const edm::EventSetup es,
const HitDoublets doublets,
const RecHitsSortedInPhi **  thirdHitMap,
const std::vector< const DetLayer * > &  thirdLayerDetLayer,
const int  nThirdLayers 
)
overridevirtual

Implements HitTripletGeneratorFromPairAndLayers.

Definition at line 174 of file PixelTripletNoTipGenerator.cc.

References Exception.

182 {
183  throw cms::Exception("Error")<<"PixelTripletNoTipGenerator::hitTriplets is not implemented \n";
184 }

Member Data Documentation

float PixelTripletNoTipGenerator::extraHitPhiToleranceForPreFiltering
private

Definition at line 36 of file PixelTripletNoTipGenerator.h.

Referenced by hitTriplets().

float PixelTripletNoTipGenerator::extraHitRPhitolerance
private

Definition at line 35 of file PixelTripletNoTipGenerator.h.

float PixelTripletNoTipGenerator::extraHitRZtolerance
private

Definition at line 34 of file PixelTripletNoTipGenerator.h.

Referenced by hitTriplets().

double PixelTripletNoTipGenerator::theChi2Cut
private

Definition at line 38 of file PixelTripletNoTipGenerator.h.

Referenced by hitTriplets().

double PixelTripletNoTipGenerator::theNSigma
private

Definition at line 37 of file PixelTripletNoTipGenerator.h.

Referenced by hitTriplets().