CMS 3D CMS Logo

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

#include <InnerDeltaPhi.h>

Public Types

typedef Basic2DVector< float > Point2D
 

Public Member Functions

 InnerDeltaPhi (const DetLayer &outlayer, const DetLayer &layer, const TrackingRegion &region, const edm::EventSetup &iSetup, bool precise=true, float extraTolerance=0.f)
 
PixelRecoRange< float > operator() (float xHit, float yHit, float zHit, float errRPhi) const
 
bool prefilter (float xHit, float yHit) const
 

Private Member Functions

void initBarrelLayer (const DetLayer &layer)
 
void initBarrelMS (const DetLayer &outLayer)
 
void initForwardLayer (const DetLayer &layer, float zMinOrigin, float zMaxOrigin)
 
void initForwardMS (const DetLayer &outLayer)
 
PixelRecoRange< float > phiRange (const Point2D &hitXY, float zHit, float errRPhi) const
 

Private Attributes

bool innerIsBarrel
 
int ol
 
bool outerIsBarrel
 
MultipleScatteringParametrisation sigma
 
float theA
 
float theB
 
float theDeltaScatt
 
float theExtraTolerance
 
bool thePrecise
 
float thePtMin
 
float theRCurvature
 
float theRLayer
 
float theROrigin
 
float theScatt0
 
float theThickness
 
Point2D theVtx
 
float theVtxZ
 

Detailed Description

Definition at line 20 of file InnerDeltaPhi.h.

Member Typedef Documentation

Definition at line 23 of file InnerDeltaPhi.h.

Constructor & Destructor Documentation

InnerDeltaPhi::InnerDeltaPhi ( const DetLayer outlayer,
const DetLayer layer,
const TrackingRegion region,
const edm::EventSetup iSetup,
bool  precise = true,
float  extraTolerance = 0.f 
)

Definition at line 82 of file InnerDeltaPhi.cc.

References PixelRecoUtilities::bendingRadius(), initBarrelLayer(), initBarrelMS(), initForwardLayer(), initForwardMS(), innerIsBarrel, TrackingRegion::originZBound(), outerIsBarrel, thePtMin, theRCurvature, and theVtxZ.

85  :
86  innerIsBarrel(layer.isBarrel()),
87  outerIsBarrel(outlayer.isBarrel()),
89  ol( outlayer.seqNum()),
90  theROrigin(region.originRBound()),
91  theRLayer(0),
92  theThickness(0),
93  theExtraTolerance(extraTolerance),
94  theA(0),
95  theB(0),
96  theVtxZ(region.origin().z()),
97  thePtMin(region.ptMin()),
98  theVtx(region.origin().x(),region.origin().y()),
99  sigma(&layer,iSetup)
100 
101 {
102  float zMinOrigin = theVtxZ-region.originZBound();
103  float zMaxOrigin = theVtxZ+region.originZBound();
105 
106 
107  if (innerIsBarrel) initBarrelLayer( layer);
108  else initForwardLayer( layer, zMinOrigin, zMaxOrigin);
109 
110  if(outerIsBarrel) initBarrelMS(outlayer);
111  else initForwardMS(outlayer);
112 
113 }
float originRBound() const
bounds the particle vertex in the transverse plane
GlobalPoint const & origin() const
void initBarrelMS(const DetLayer &outLayer)
T y() const
Definition: PV3DBase.h:63
Point2D theVtx
Definition: InnerDeltaPhi.h:61
int seqNum() const
Definition: DetLayer.h:36
T z() const
Definition: PV3DBase.h:64
float originZBound() const
bounds the particle vertex in the longitudinal plane
bool isBarrel() const
Definition: DetLayer.h:32
void initForwardMS(const DetLayer &outLayer)
float theThickness
Definition: InnerDeltaPhi.h:49
void initBarrelLayer(const DetLayer &layer)
float theRCurvature
Definition: InnerDeltaPhi.h:53
void initForwardLayer(const DetLayer &layer, float zMinOrigin, float zMaxOrigin)
float ptMin() const
minimal pt of interest
T bendingRadius(T pt, const edm::EventSetup &iSetup)
float theExtraTolerance
Definition: InnerDeltaPhi.h:54
T x() const
Definition: PV3DBase.h:62
MultipleScatteringParametrisation sigma
Definition: InnerDeltaPhi.h:64

Member Function Documentation

void InnerDeltaPhi::initBarrelLayer ( const DetLayer layer)
private

Definition at line 146 of file InnerDeltaPhi.cc.

References Surface::bounds(), BarrelDetLayer::specificSurface(), GeometricSearchDet::surface(), theRLayer, theThickness, and Bounds::thickness().

Referenced by InnerDeltaPhi().

147 {
148  const BarrelDetLayer& bl = static_cast<const BarrelDetLayer&>(layer);
149  float rLayer = bl.specificSurface().radius();
150 
151  // the maximal delta phi will be for the innermost hits
152  theThickness = layer.surface().bounds().thickness();
153  theRLayer = rLayer - 0.5f*theThickness;
154 }
const Bounds & bounds() const
Definition: Surface.h:120
virtual const BoundCylinder & specificSurface() const final
Extension of the interface.
virtual const BoundSurface & surface() const =0
The surface of the GeometricSearchDet.
float theThickness
Definition: InnerDeltaPhi.h:49
virtual float thickness() const =0
void InnerDeltaPhi::initBarrelMS ( const DetLayer outLayer)
private

Definition at line 116 of file InnerDeltaPhi.cc.

References Surface::bounds(), Bounds::length(), ol, sigma, BarrelDetLayer::specificSurface(), GeometricSearchDet::surface(), theDeltaScatt, thePtMin, and theScatt0.

Referenced by InnerDeltaPhi().

116  {
117  const BarrelDetLayer& bl = static_cast<const BarrelDetLayer&>(outLayer);
118  float rLayer = bl.specificSurface().radius();
119  auto zmax = 0.5f*outLayer.surface().bounds().length();
120  PixelRecoPointRZ zero(0., 0.);
121  PixelRecoPointRZ point1(rLayer, 0.);
122  PixelRecoPointRZ point2(rLayer, zmax);
123  auto scatt1 = 3.f*sigma(thePtMin,zero, point1, ol);
124  auto scatt2 = 3.f*sigma(thePtMin,zero, point2, ol);
125  theDeltaScatt = (scatt2-scatt1)/zmax;
126  theScatt0 = scatt1;
127 }
virtual float length() const =0
const Bounds & bounds() const
Definition: Surface.h:120
virtual const BoundCylinder & specificSurface() const final
Extension of the interface.
float theDeltaScatt
Definition: InnerDeltaPhi.h:51
virtual const BoundSurface & surface() const =0
The surface of the GeometricSearchDet.
MultipleScatteringParametrisation sigma
Definition: InnerDeltaPhi.h:64
void InnerDeltaPhi::initForwardLayer ( const DetLayer layer,
float  zMinOrigin,
float  zMaxOrigin 
)
private

Definition at line 156 of file InnerDeltaPhi.cc.

References Surface::bounds(), GeometricSearchDet::position(), ForwardDetLayer::specificSurface(), GeometricSearchDet::surface(), theA, theB, theRLayer, theThickness, Bounds::thickness(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by InnerDeltaPhi().

158 {
159  const ForwardDetLayer &fl = static_cast<const ForwardDetLayer&>(layer);
160  theRLayer = fl.specificSurface().innerRadius();
161  float layerZ = layer.position().z();
162  theThickness = layer.surface().bounds().thickness();
163  float layerZmin = layerZ > 0 ? layerZ-0.5f*theThickness: layerZ+0.5f*theThickness;
164  theB = layerZ > 0 ? zMaxOrigin : zMinOrigin;
165  theA = layerZmin - theB;
166 }
const Bounds & bounds() const
Definition: Surface.h:120
T z() const
Definition: PV3DBase.h:64
virtual const BoundSurface & surface() const =0
The surface of the GeometricSearchDet.
virtual const BoundDisk & specificSurface() const final
float theThickness
Definition: InnerDeltaPhi.h:49
virtual float thickness() const =0
virtual const Surface::PositionType & position() const
Returns position of the surface.
void InnerDeltaPhi::initForwardMS ( const DetLayer outLayer)
private

Definition at line 129 of file InnerDeltaPhi.cc.

References ol, GeometricSearchDet::position(), sigma, ForwardDetLayer::specificSurface(), theDeltaScatt, thePtMin, theScatt0, theVtxZ, and PV3DBase< T, PVType, FrameType >::z().

Referenced by InnerDeltaPhi().

129  {
130  const ForwardDetLayer &fl = static_cast<const ForwardDetLayer&>(outLayer);
131  auto minR = fl.specificSurface().innerRadius();
132  auto maxR = fl.specificSurface().outerRadius();
133  auto layerZ = outLayer.position().z();
134  // compute min and max multiple scattering correction
135  PixelRecoPointRZ zero(0., theVtxZ);
136  PixelRecoPointRZ point1(minR, layerZ);
137  PixelRecoPointRZ point2(maxR, layerZ);
138  auto scatt1 = 3.f*sigma(thePtMin,zero, point1, ol);
139  auto scatt2 = 3.f*sigma(thePtMin,zero, point2, ol);
140  theDeltaScatt = (scatt2-scatt1)/(maxR-minR);
141  theScatt0 = scatt1 - theDeltaScatt*minR;
142 }
float theDeltaScatt
Definition: InnerDeltaPhi.h:51
T z() const
Definition: PV3DBase.h:64
virtual const BoundDisk & specificSurface() const final
virtual const Surface::PositionType & position() const
Returns position of the surface.
MultipleScatteringParametrisation sigma
Definition: InnerDeltaPhi.h:64
PixelRecoRange<float> InnerDeltaPhi::operator() ( float  xHit,
float  yHit,
float  zHit,
float  errRPhi 
) const
inline

Definition at line 36 of file InnerDeltaPhi.h.

36  {
37  return phiRange( Point2D(xHit,yHit), zHit, errRPhi);
38  }
PixelRecoRange< float > phiRange(const Point2D &hitXY, float zHit, float errRPhi) const
Basic2DVector< float > Point2D
Definition: InnerDeltaPhi.h:23
PixelRecoRange< float > InnerDeltaPhi::phiRange ( const Point2D hitXY,
float  zHit,
float  errRPhi 
) const
private

Definition at line 170 of file InnerDeltaPhi.cc.

References funct::abs(), hiPixelPairStep_cff::deltaPhi, Basic2DVector< T >::dot(), dt, MillePedeFileConverter_cfg::e, f, objects.autophobj::float, innerIsBarrel, M_PI, Basic2DVector< T >::mag(), margin, SiStripPI::max, min(), pileupDistInMC::num, outerIsBarrel, funct::sqr(), mathSSE::sqrt(), lumiQTWidget::t, theA, theB, theDeltaScatt, theExtraTolerance, thePrecise, theRCurvature, theRLayer, theROrigin, theScatt0, theThickness, theVtx, Basic2DVector< T >::unit(), and w.

171 {
172  float rLayer = theRLayer;
173  Point2D crossing;
174 
175  Point2D dHit = hitXY-theVtx;
176  auto dHitmag = dHit.mag();
177  float dLayer = 0.;
178  float dL = 0.;
179 
180  // track is crossing layer with angle such as:
181  // this factor should be taken in computation of eror projection
182  float cosCross = 0;
183 
184 
185  //
186  // compute crossing of stright track with inner layer
187  //
188  if (!innerIsBarrel) {
189  auto t = theA/(hitZ-theB); auto dt = std::abs(theThickness/(hitZ-theB));
190  crossing = theVtx + t*dHit;
191  rLayer = crossing.mag();
192  dLayer = t*dHitmag; dL = dt * dHitmag;
193  cosCross = std::abs( dHit.unit().dot(crossing.unit()));
194  } else {
195 
196  //
197  // compute crossing of track with layer
198  // dHit - from VTX to outer hit
199  // rLayer - layer radius
200  // dLayer - distance from VTX to inner layer in direction of dHit
201  // vect(rLayer) = vect(rVTX) + vect(dHit).unit * dLayer
202  // rLayer^2 = (vect(rVTX) + vect(dHit).unit * dLayer)^2 and we have square eqation for dLayer
203  //
204  // barrel case
205  //
206  auto vtxmag2 = theVtx.mag2();
207  if (vtxmag2 < 1.e-10f) {
208  dLayer = rLayer;
209  }
210  else {
211  // there are cancellation here....
212  double var_c = vtxmag2-sqr(rLayer);
213  double var_b = theVtx.dot(dHit.unit());
214  double var_delta = sqr(var_b)-var_c;
215  if (var_delta <=0.) var_delta = 0;
216  dLayer = -var_b + std::sqrt(var_delta); //only the value along vector is OK.
217  }
218  crossing = theVtx+ dHit.unit() * dLayer;
219  cosCross = std::abs( dHit.unit().dot(crossing.unit()));
220  dL = theThickness/cosCross;
221  }
222 
223 #ifdef USE_VECTORS_HERE
224  cms_float32x4_t num{dHitmag,dLayer,theROrigin * (dHitmag-dLayer),1.f};
225  cms_float32x4_t den{2*theRCurvature,2*theRCurvature,dHitmag*dLayer,1.f};
226  auto phis = f_asin07f(num/den);
227  phis = phis*dLayer/(rLayer*cosCross);
228  auto deltaPhi = std::abs(phis[0]-phis[1]);
229  auto deltaPhiOrig = phis[2];
230 #else
231 #warning no vector!
232  auto alphaHit = cropped_asin( dHitmag/(2*theRCurvature));
233  auto OdeltaPhi = std::abs( alphaHit - cropped_asin( dLayer/(2*theRCurvature)));
234  OdeltaPhi *= dLayer/(rLayer*cosCross);
235  // compute additional delta phi due to origin radius
236  auto OdeltaPhiOrig = cropped_asin( theROrigin * (dHitmag-dLayer) / (dHitmag*dLayer));
237  OdeltaPhiOrig *= dLayer/(rLayer*cosCross);
238  // std::cout << "dphi " << OdeltaPhi<<'/'<<OdeltaPhiOrig << ' ' << deltaPhi<<'/'<<deltaPhiOrig << std::endl;
239 
240  auto deltaPhi = OdeltaPhi;
241  auto deltaPhiOrig = OdeltaPhiOrig;
242 #endif
243 
244  // additinal angle due to not perpendicular stright line crossing (for displaced beam)
245  // double dPhiCrossing = (cosCross > 0.9999) ? 0 : dL * sqrt(1-sqr(cosCross))/ rLayer;
246  Point2D crossing2 = theVtx + dHit.unit()* (dLayer+dL);
247  auto phicross2 = f_phi(crossing2);
248  auto phicross1 = f_phi(crossing);
249  auto dphicross = phicross2-phicross1;
250  if (dphicross < -float(M_PI)) dphicross += float(2*M_PI);
251  if (dphicross > float(M_PI)) dphicross -= float(2*M_PI);
252  if (dphicross > float(M_PI/2)) dphicross = 0.; // something wrong?
253  phicross2 = phicross1 + dphicross;
254 
255 
256 
257  // inner hit error taken as constant
258  auto deltaPhiHit = theExtraTolerance / rLayer;
259 
260  // outer hit error
261  // double deltaPhiHitOuter = errRPhi/rLayer;
262  auto deltaPhiHitOuter = errRPhi/hitXY.mag();
263 
264  auto margin = deltaPhi+deltaPhiOrig+deltaPhiHit+deltaPhiHitOuter ;
265 
266  if (thePrecise) {
267  // add multiple scattering correction
268 
269  /*
270  PixelRecoPointRZ zero(0., theVtxZ);
271  PixelRecoPointRZ point(hitXY.mag(), hitZ);
272  auto scatt = 3.f*sigma(thePtMin,zero, point, ol);
273  */
274 
275  auto w = outerIsBarrel ? std::abs(hitZ) : hitXY.mag();
276  auto nscatt = theScatt0 + theDeltaScatt*w;
277 
278  // std::cout << "scatt " << (outerIsBarrel ? "B" : "F") << (innerIsBarrel ? "B " : "F ")
279  // << scatt << ' ' << nscatt << ' ' << nscatt/scatt << std::endl;
280 
281  margin += nscatt/ rLayer ;
282  }
283 
284  return PixelRecoRange<float>( std::min(phicross1,phicross2)-margin,
285  std::max(phicross1,phicross2)+margin);
286 }
float dt
Definition: AMPTWrapper.h:126
const double w
Definition: UKUtility.cc:23
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Point2D theVtx
Definition: InnerDeltaPhi.h:61
Basic2DVector unit() const
float theDeltaScatt
Definition: InnerDeltaPhi.h:51
T sqrt(T t)
Definition: SSEVec.h:18
Double_t margin
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
T min(T a, T b)
Definition: MathUtil.h:58
#define M_PI
float theThickness
Definition: InnerDeltaPhi.h:49
float theRCurvature
Definition: InnerDeltaPhi.h:53
Square< F >::type sqr(const F &f)
Definition: Square.h:13
float theExtraTolerance
Definition: InnerDeltaPhi.h:54
bool InnerDeltaPhi::prefilter ( float  xHit,
float  yHit 
) const
inline

Definition at line 32 of file InnerDeltaPhi.h.

Referenced by HitPairGeneratorFromLayerPair::doublets().

32  {
33  return xHit*xHit + yHit*yHit > theRLayer*theRLayer;
34  }

Member Data Documentation

bool InnerDeltaPhi::innerIsBarrel
private

Definition at line 42 of file InnerDeltaPhi.h.

Referenced by InnerDeltaPhi(), and phiRange().

int InnerDeltaPhi::ol
private

Definition at line 45 of file InnerDeltaPhi.h.

Referenced by initBarrelMS(), and initForwardMS().

bool InnerDeltaPhi::outerIsBarrel
private

Definition at line 43 of file InnerDeltaPhi.h.

Referenced by InnerDeltaPhi(), and phiRange().

MultipleScatteringParametrisation InnerDeltaPhi::sigma
private

Definition at line 64 of file InnerDeltaPhi.h.

Referenced by initBarrelMS(), and initForwardMS().

float InnerDeltaPhi::theA
private

Definition at line 55 of file InnerDeltaPhi.h.

Referenced by initForwardLayer(), and phiRange().

float InnerDeltaPhi::theB
private

Definition at line 56 of file InnerDeltaPhi.h.

Referenced by initForwardLayer(), and phiRange().

float InnerDeltaPhi::theDeltaScatt
private

Definition at line 51 of file InnerDeltaPhi.h.

Referenced by initBarrelMS(), initForwardMS(), and phiRange().

float InnerDeltaPhi::theExtraTolerance
private

Definition at line 54 of file InnerDeltaPhi.h.

Referenced by phiRange().

bool InnerDeltaPhi::thePrecise
private

Definition at line 44 of file InnerDeltaPhi.h.

Referenced by phiRange().

float InnerDeltaPhi::thePtMin
private

Definition at line 59 of file InnerDeltaPhi.h.

Referenced by initBarrelMS(), initForwardMS(), and InnerDeltaPhi().

float InnerDeltaPhi::theRCurvature
private

Definition at line 53 of file InnerDeltaPhi.h.

Referenced by InnerDeltaPhi(), and phiRange().

float InnerDeltaPhi::theRLayer
private

Definition at line 48 of file InnerDeltaPhi.h.

Referenced by initBarrelLayer(), initForwardLayer(), and phiRange().

float InnerDeltaPhi::theROrigin
private

Definition at line 47 of file InnerDeltaPhi.h.

Referenced by phiRange().

float InnerDeltaPhi::theScatt0
private

Definition at line 50 of file InnerDeltaPhi.h.

Referenced by initBarrelMS(), initForwardMS(), and phiRange().

float InnerDeltaPhi::theThickness
private

Definition at line 49 of file InnerDeltaPhi.h.

Referenced by initBarrelLayer(), initForwardLayer(), and phiRange().

Point2D InnerDeltaPhi::theVtx
private

Definition at line 61 of file InnerDeltaPhi.h.

Referenced by phiRange().

float InnerDeltaPhi::theVtxZ
private

Definition at line 58 of file InnerDeltaPhi.h.

Referenced by initForwardMS(), and InnerDeltaPhi().