CMS 3D CMS Logo

Public Types | Public Member Functions | Private Member Functions | Private Attributes

InnerDeltaPhi Class Reference

#include <InnerDeltaPhi.h>

List of all members.

Public Types

typedef Basic2DVector< double > Point2D

Public Member Functions

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

Private Member Functions

void initBarrelLayer (const DetLayer &layer)
void initForwardLayer (const DetLayer &layer, float zMinOrigin, float zMaxOrigin)
float minRadius (float hitR, float hitZ) const
PixelRecoRange< float > phiRange (const Point2D &hitXY, float zHit, float errRPhi) const

Private Attributes

MultipleScatteringParametrisation sigma
float theA
float theB
float theExtraTolerance
bool thePrecise
float thePtMin
float theRCurvature
bool theRDefined
float theRLayer
float theROrigin
float theThickness
Point2D theVtx
float theVtxZ

Detailed Description

Definition at line 16 of file InnerDeltaPhi.h.


Member Typedef Documentation

Definition at line 19 of file InnerDeltaPhi.h.


Constructor & Destructor Documentation

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

Definition at line 30 of file InnerDeltaPhi.cc.

References Reference_intrackfit_cff::barrel, PixelRecoUtilities::bendingRadius(), initBarrelLayer(), initForwardLayer(), DetLayer::location(), TrackingRegion::originZBound(), thePtMin, theRCurvature, and theVtxZ.

  : thePrecise(precise),
    theROrigin(region.originRBound()),
    theRLayer(0),
    theThickness(0),
    theExtraTolerance(extraTolerance),
    theA(0),
    theB(0),
    theVtxZ(region.origin().z()),
    thePtMin(region.ptMin()),
    theVtx(region.origin().x(),region.origin().y()),
    sigma(&layer,iSetup)
{
  float zMinOrigin = theVtxZ-region.originZBound();
  float zMaxOrigin = theVtxZ+region.originZBound();
  theRCurvature = PixelRecoUtilities::bendingRadius(thePtMin,iSetup);
 

  if (layer.location() == GeomDetEnumerators::barrel) initBarrelLayer( layer);
  else initForwardLayer( layer, zMinOrigin, zMaxOrigin);

}
InnerDeltaPhi::~InnerDeltaPhi ( )

Definition at line 57 of file InnerDeltaPhi.cc.

{}

Member Function Documentation

void InnerDeltaPhi::initBarrelLayer ( const DetLayer layer) [private]

Definition at line 59 of file InnerDeltaPhi.cc.

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

Referenced by InnerDeltaPhi().

{
  const BarrelDetLayer& bl = static_cast<const BarrelDetLayer&>(layer); 
  float rLayer = bl.specificSurface().radius(); 

  // the maximal delta phi will be for the innermost hits
  theThickness = layer.surface().bounds().thickness();
  theRLayer = rLayer - 0.5f*theThickness;
  theRDefined = true;
}
void InnerDeltaPhi::initForwardLayer ( const DetLayer layer,
float  zMinOrigin,
float  zMaxOrigin 
) [private]

Definition at line 70 of file InnerDeltaPhi.cc.

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

Referenced by InnerDeltaPhi().

{
  const ForwardDetLayer &fl = static_cast<const ForwardDetLayer&>(layer);
  theRLayer = fl.specificSurface().innerRadius();
  float layerZ = layer.position().z();
  theThickness = layer.surface().bounds().thickness();
  float layerZmin = layerZ > 0 ? layerZ-0.5f*theThickness: layerZ+0.5f*theThickness;
  theB = layerZ > 0 ? zMaxOrigin : zMinOrigin;
  theA = layerZmin - theB;
  theRDefined = false;
}
float InnerDeltaPhi::minRadius ( float  hitR,
float  hitZ 
) const [private]

Definition at line 229 of file InnerDeltaPhi.cc.

References max(), theA, theB, theRDefined, and theRLayer.

Referenced by operator()().

{
  if (theRDefined) return theRLayer;
  else {
    float rmin = (theA*hitR)/(hitZ-theB);
    return ( rmin> 0) ? std::max( rmin, theRLayer) : theRLayer;
  }
}
PixelRecoRange< float > InnerDeltaPhi::operator() ( float  rHit,
float  phiHit,
float  zHit,
float  errRPhi 
) const

Definition at line 219 of file InnerDeltaPhi.cc.

References funct::cos(), phiRange(), and funct::sin().

{
//     float phiM =  operator()( rHit,zHit,errRPhi);
//     return PixelRecoRange<float>(phiHit-phiM,phiHit+phiM);

  Point2D hitXY( rHit*std::cos(phiHit), rHit*std::sin(phiHit));
  return phiRange(hitXY,zHit,errRPhi);
}
float InnerDeltaPhi::operator() ( float  rHit,
float  zHit,
float  errRPhi 
) const

Definition at line 186 of file InnerDeltaPhi.cc.

References abs, SiPixelRawToDigiRegional_cfi::deltaPhi, minRadius(), point, sigma, theExtraTolerance, thePrecise, thePtMin, theRCurvature, theROrigin, theVtxZ, and zero.

{
  // alpha - angle between particle direction at vertex and position of hit.
  // (pi/2 - alpha) - angle hit-vertex-cernter_of_curvature
  // cos (pi/2 - alpha) = (hRhi/2) / theRCurvature
  // so:

  float alphaHit = std::asin( rHit/(2*theRCurvature));

  
  float rMin = minRadius( rHit, zHit);
  float deltaPhi = std::abs( alphaHit - std::asin( rMin/(2*theRCurvature)));

  // compute additional delta phi due to origin radius
  float deltaPhiOrig = asin( theROrigin * (rHit-rMin) / (rHit*rMin));

  // hit error taken as constant
  float deltaPhiHit = theExtraTolerance / rMin;

  if (!thePrecise) {
    return deltaPhi+deltaPhiOrig+deltaPhiHit;
  } else {
    // add multiple scattering correction
    PixelRecoPointRZ zero(0., theVtxZ);
    PixelRecoPointRZ point(rHit, zHit);
    float scatt = 3.f*sigma(thePtMin,zero, point) / rMin; 
    float deltaPhiHitOuter = errRPhi/rMin; 
   
    return deltaPhi+deltaPhiOrig+deltaPhiHit + scatt + deltaPhiHitOuter;
  }

}
PixelRecoRange< float > InnerDeltaPhi::phiRange ( const Point2D hitXY,
float  zHit,
float  errRPhi 
) const [private]

Definition at line 85 of file InnerDeltaPhi.cc.

References abs, cropped_asin(), SiPixelRawToDigiRegional_cfi::deltaPhi, Basic2DVector< T >::dot(), dt, alignCSCRings::e, M_PI, Basic2DVector< T >::mag(), max(), min, Basic2DVector< T >::phi(), point, sigma, funct::sqr(), mathSSE::sqrt(), lumiQTWidget::t, theA, theB, theExtraTolerance, thePrecise, thePtMin, theRCurvature, theRDefined, theRLayer, theROrigin, theThickness, theVtx, theVtxZ, Basic2DVector< T >::unit(), and zero.

Referenced by operator()().

{
  double rLayer = theRLayer;
  bool checkCrossing = true;
  Point2D crossing;

  Point2D dHit = hitXY-theVtx;
  double  dHitmag = dHit.mag();
  double  dLayer = 0.;
  double dL = 0.;
  //
  // compute crossing of stright track with inner layer
  //
  if (!theRDefined) {
    double t = theA/(hitZ-theB); double dt = std::abs(theThickness/(hitZ-theB));
    crossing = theVtx + t*dHit;
    rLayer =  crossing.mag();
    dLayer = t*dHitmag;           dL = dt * dHitmag; 
    checkCrossing = false;
    if (rLayer < theRLayer) {
      checkCrossing = true;
      rLayer = theRLayer;
      dL = 0.;
    } 
  }

  //
  // compute crossing of track with layer
  // dHit - from VTX to outer hit
  // rLayer - layer radius
  // dLayer - distance from VTX to inner layer in direction of dHit
  // vect(rLayer) = vect(rVTX) + vect(dHit).unit * dLayer
  //     rLayer^2 = (ect(rVTX) + vect(dHit).unit * dLayer)^2 and we have square eqation for dLayer 
  //
  // barrel case
  //

  if (checkCrossing) {
    double vtxmag2 = theVtx.mag2();
    if (vtxmag2 < 1.e-10) {
      dLayer = rLayer;
    }
    else { 
      double var_c = theVtx.mag2()-sqr(rLayer);
      double var_b = 2*theVtx.dot(dHit.unit());
      double var_delta = sqr(var_b)-4*var_c;
      if (var_delta <=0.) var_delta = 0;
      dLayer = 0.5*(-var_b + std::sqrt(var_delta)); //only the value along vector is OK. 
    }
    crossing = theVtx+ dHit.unit() * dLayer;
    double cosCross = std::abs( dHit.unit().dot(crossing.unit()));
    dL = theThickness/cosCross; 
  }


  // track is crossing layer with angle such as:
  // this factor should be taken in computation of eror projection
  double cosCross = std::abs( dHit.unit().dot(crossing.unit()));

  double alphaHit = cropped_asin( dHitmag/(2*theRCurvature)); 
  double deltaPhi = fabs( alphaHit - cropped_asin( dLayer/(2*theRCurvature)));
  deltaPhi *= dLayer/(rLayer*cosCross);  

  // additinal angle due to not perpendicular stright line crossing  (for displaced beam)
  //  double dPhiCrossing = (cosCross > 0.9999) ? 0 : dL *  sqrt(1-sqr(cosCross))/ rLayer;
  Point2D crossing2 = theVtx + dHit.unit()* (dLayer+dL);
  double phicross2 = crossing2.phi();  
  double phicross1 = crossing.phi();
  double dphicross = phicross2-phicross1;
  if (dphicross < -M_PI) dphicross += 2*M_PI;
  if (dphicross >  M_PI) dphicross -= 2*M_PI;
  if (dphicross > M_PI/2) dphicross = 0.;  // something wrong?
  phicross2 = phicross1 + dphicross;
        

  // compute additional delta phi due to origin radius
  double deltaPhiOrig = cropped_asin( theROrigin * (dHitmag-dLayer) / (dHitmag*dLayer));
  deltaPhiOrig *= dLayer/(rLayer*cosCross);

  // inner hit error taken as constant
  double deltaPhiHit = theExtraTolerance / rLayer;

  // outer hit error
//   double deltaPhiHitOuter = errRPhi/rLayer; 
    double deltaPhiHitOuter = errRPhi/hitXY.mag();

  double margin = deltaPhi+deltaPhiOrig+deltaPhiHit+deltaPhiHitOuter ;

  if (thePrecise) {
    // add multiple scattering correction
    PixelRecoPointRZ zero(0., theVtxZ);
    PixelRecoPointRZ point(hitXY.mag(), hitZ);
    double scatt = 3.f*sigma(thePtMin,zero, point) / rLayer; 
   
    margin += scatt ;
  }
  
  return PixelRecoRange<float>( std::min(phicross1,phicross2)-margin, 
                                std::max(phicross1,phicross2)+margin);
}

Member Data Documentation

Definition at line 54 of file InnerDeltaPhi.h.

Referenced by operator()(), and phiRange().

float InnerDeltaPhi::theA [private]

Definition at line 45 of file InnerDeltaPhi.h.

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

float InnerDeltaPhi::theB [private]

Definition at line 46 of file InnerDeltaPhi.h.

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

Definition at line 44 of file InnerDeltaPhi.h.

Referenced by operator()(), and phiRange().

bool InnerDeltaPhi::thePrecise [private]

Definition at line 36 of file InnerDeltaPhi.h.

Referenced by operator()(), and phiRange().

float InnerDeltaPhi::thePtMin [private]

Definition at line 49 of file InnerDeltaPhi.h.

Referenced by InnerDeltaPhi(), operator()(), and phiRange().

Definition at line 43 of file InnerDeltaPhi.h.

Referenced by InnerDeltaPhi(), operator()(), and phiRange().

Definition at line 35 of file InnerDeltaPhi.h.

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

float InnerDeltaPhi::theRLayer [private]

Definition at line 40 of file InnerDeltaPhi.h.

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

float InnerDeltaPhi::theROrigin [private]

Definition at line 39 of file InnerDeltaPhi.h.

Referenced by operator()(), and phiRange().

float InnerDeltaPhi::theThickness [private]

Definition at line 41 of file InnerDeltaPhi.h.

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

Definition at line 51 of file InnerDeltaPhi.h.

Referenced by phiRange().

float InnerDeltaPhi::theVtxZ [private]

Definition at line 48 of file InnerDeltaPhi.h.

Referenced by InnerDeltaPhi(), operator()(), and phiRange().