CMS 3D CMS Logo

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

#include <TrackFitter.h>

Inheritance diagram for TrackFitter:
PixelFitterBase

Public Member Functions

std::unique_ptr< reco::Trackrun (const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region) const override
 
 TrackFitter (const edm::EventSetup *es, const TrackerGeometry *tracker, const MagneticField *field, const TransientTrackingRecHitBuilder *ttrhBuilder)
 
virtual ~TrackFitter ()
 
- Public Member Functions inherited from PixelFitterBase
virtual reco::Trackrun (const edm::EventSetup &es, const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region) const
 
virtual reco::Trackrun (const edm::Event &ev, const edm::EventSetup &es, const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region) const
 
virtual ~PixelFitterBase ()
 

Private Member Functions

float getCotThetaAndUpdateZip (const GlobalPoint &inner, const GlobalPoint &outer, float radius, float phi, float d0, float &zip) const
 
void getErrTipAndErrZip (float pt, float eta, float &errZip, float &errTip) const
 
float getPhi (float xC, float yC, int charge) const
 
float getZip (float d0, float curv, const GlobalPoint &inner, const GlobalPoint &outer) const
 

Private Attributes

const edm::EventSetuptheES
 
const MagneticFieldtheField
 
const TrackerGeometrytheTracker
 
const TransientTrackingRecHitBuildertheTTRecHitBuilder
 

Detailed Description

Definition at line 17 of file TrackFitter.h.

Constructor & Destructor Documentation

TrackFitter::TrackFitter ( const edm::EventSetup es,
const TrackerGeometry tracker,
const MagneticField field,
const TransientTrackingRecHitBuilder ttrhBuilder 
)
inline

Definition at line 20 of file TrackFitter.h.

21  :
22  theES(es), theTracker(tracker), theField(field), theTTRecHitBuilder(ttrhBuilder)
23  {}
const TrackerGeometry * theTracker
Definition: TrackFitter.h:39
const edm::EventSetup * theES
Definition: TrackFitter.h:38
const TransientTrackingRecHitBuilder * theTTRecHitBuilder
Definition: TrackFitter.h:41
const MagneticField * theField
Definition: TrackFitter.h:40
virtual TrackFitter::~TrackFitter ( )
inlinevirtual

Member Function Documentation

float TrackFitter::getCotThetaAndUpdateZip ( const GlobalPoint inner,
const GlobalPoint outer,
float  radius,
float  phi,
float  d0,
float &  zip 
) const
private

Definition at line 130 of file TrackFitter.cc.

References funct::cos(), runTauDisplay::dr, allConversions_cfi::dz, getPhi(), listHistos::IP, M_PI_2, perp(), funct::sin(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by run(), and ~TrackFitter().

132 {
133  float chi = phi - M_PI_2;
134  GlobalPoint IP(d0*cos(chi), d0*sin(chi),zip);
135 
136  float phi1 = 2*asin(0.5*(inner - IP).perp()/radius);
137  float phi2 = 2*asin(0.5*(outer - IP).perp()/radius);
138 
139  float dr = radius*(phi2 - phi1);
140  float dz = outer.z()-inner.z();
141 
142  // Recalculate ZIP
143  zip = (inner.z()*phi2 - outer.z()*phi1)/(phi2 - phi1);
144 
145  return (fabs(dr) > 1.e-3) ? dz/dr : 0;
146 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define M_PI_2
OutputIterator zip(InputIterator1 first1, InputIterator1 last1, InputIterator2 first2, InputIterator2 last2, OutputIterator result, Compare comp)
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
T perp() const
Magnitude of transverse component.
void TrackFitter::getErrTipAndErrZip ( float  pt,
float  eta,
float &  errZip,
float &  errTip 
) const
private

Definition at line 182 of file TrackFitter.cc.

References EnergyCorrector::pt, and mathSSE::sqrt().

Referenced by getZip(), and ~TrackFitter().

183 {
184  float coshEta = cosh(eta);
185 
186  { // transverse
187  float c_ms = 0.0115; //0.0115;
188  float s_le = 0.0095; //0.0123;
189  float s_ms2 = c_ms*c_ms / (pt*pt) * coshEta;
190 
191  errTip = sqrt(s_le*s_le + s_ms2 );
192  }
193 
194  { // z
195  float c_ms = 0.0070;
196  float s_le = 0.0135;
197 
198  errZip = sqrt( (s_le*s_le + c_ms*c_ms/(pt*pt)) * coshEta*coshEta*coshEta);
199  }
200 }
T sqrt(T t)
Definition: SSEVec.h:18
float TrackFitter::getPhi ( float  xC,
float  yC,
int  charge 
) const
private

Definition at line 150 of file TrackFitter.cc.

References getZip().

Referenced by getCotThetaAndUpdateZip(), and ~TrackFitter().

151 {
152  float phiC;
153 
154  if (charge>0) phiC = atan2(xC,-yC);
155  else phiC = atan2(-xC,yC);
156 
157  return phiC;
158 }
float TrackFitter::getZip ( float  d0,
float  curv,
const GlobalPoint inner,
const GlobalPoint outer 
) const
private

Definition at line 162 of file TrackFitter.cc.

References getErrTipAndErrZip(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::perp2(), diffTwoXMLs::r1, diffTwoXMLs::r2, and PV3DBase< T, PVType, FrameType >::z().

Referenced by getPhi(), and ~TrackFitter().

164 {
165  // phi = asin(r*rho/2) with asin(x) ~= x+x**3/(2*3)
166  float rho3 = curv*curv*curv;
167 
168  float r1 = inner.perp();
169  double phi1 = r1*curv/2 + inner.perp2()*r1*rho3/48.;
170 
171  float r2 = outer.perp();
172  double phi2 = r2*curv/2 + outer.perp2()*r2*rho3/48.;
173 
174  double z1 = inner.z();
175  double z2 = outer.z();
176 
177  return z1 - phi1/(phi1-phi2)*(z1-z2);
178 }
T perp() const
Definition: PV3DBase.h:72
T perp2() const
Definition: PV3DBase.h:71
T z() const
Definition: PV3DBase.h:64
std::unique_ptr< reco::Track > TrackFitter::run ( const std::vector< const TrackingRecHit * > &  hits,
const TrackingRegion region 
) const
overridevirtual

Reimplemented from PixelFitterBase.

Definition at line 53 of file TrackFitter.cc.

References PixelTrackBuilder::build(), CircleFromThreePoints::center(), ALCARECOTkAlJpsiMuMu_cff::charge, HiEvtPlane_cfi::chi2, RZLine::chi2(), PixelRecoUtilities::curvature(), CircleFromThreePoints::curvature(), declareDynArray, vertexPlots::e4, getCharge(), getCotThetaAndUpdateZip(), mps_fire::i, PixelRecoUtilities::inversePt(), gedGsfElectrons_cfi::isBarrel, Basic2DVector< T >::mag(), nhits, EnergyCorrector::pt, rpcPointValidation_cfi::recHit, btvTracks_cfi::tip, Basic2DVector< T >::x(), Basic2DVector< T >::y(), and ComparisonHelper::zip().

Referenced by ~TrackFitter().

55 {
56  std::unique_ptr<reco::Track> ret;
57 
58  int nhits = hits.size();
59  if(nhits <2) return ret;
60 
61  declareDynArray(GlobalPoint,nhits, points);
63  declareDynArray(bool,nhits, isBarrel);
64 
65  unsigned int i=0;
66  for (auto const & ih : hits)
67  {
68  auto recHit = theTTRecHitBuilder->build(ih);
69 
70  points[i] = recHit->globalPosition();
71  errors[i] = recHit->globalPositionError();
72  isBarrel[++i] = recHit->detUnit()->type().isBarrel();
73  }
74 
75  CircleFromThreePoints circle = (nhits==2) ?
76  CircleFromThreePoints(GlobalPoint(0.,0.,0.), points[0], points[1]) :
77  CircleFromThreePoints(points[0],points[1],points[2]);
78 
79  int charge = getCharge(points);
80  float curvature = circle.curvature();
81 
82  // pt
83  float invPt = PixelRecoUtilities::inversePt(curvature, *theES);
84  float valPt = (invPt > 1.e-4) ? 1./invPt : 1.e4;
85  float errPt = 0.055*valPt + 0.017*valPt*valPt;
86 
87  CircleFromThreePoints::Vector2D center = circle.center();
88 
89  // tip
90  float valTip = charge * (center.mag()-1/curvature);
91  // zip
92  float valZip = getZip(valTip, curvature, points[0],points[1]);
93  // phi
94  float valPhi = getPhi(center.x(), center.y(), charge);
95  // cot(theta), update zip
96  float valCotTheta =
97  getCotThetaAndUpdateZip(points[0],points[1], 1/curvature,
98  valPhi,valTip,valZip);
99 
100  // errors
101  float errTip, errZip;
102  getErrTipAndErrZip(valPt, points.back().eta(), errTip,errZip);
103  float errPhi = 0.002;
104  float errCotTheta = 0.002;
105 
106  float chi2 = 0;
107  if(nhits > 2)
108  {
109  RZLine rzLine(points,errors,isBarrel);
110  chi2 = rzLine.chi2();
111  }
112 
113  // build pixel track
114  PixelTrackBuilder builder;
115 
116  Measurement1D pt (valPt, errPt);
117  Measurement1D phi (valPhi, errPhi);
118  Measurement1D cotTheta(valCotTheta, errCotTheta);
119  Measurement1D tip (valTip, errTip);
120  Measurement1D zip (valZip, errZip);
121 
122  ret.reset(builder.build(pt, phi, cotTheta, tip, zip, chi2,
123  charge, hits, theField));
124  return ret;
125 }
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
void getErrTipAndErrZip(float pt, float eta, float &errZip, float &errTip) const
Definition: TrackFitter.cc:182
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
T inversePt(T curvature, const edm::EventSetup &iSetup)
float getZip(float d0, float curv, const GlobalPoint &inner, const GlobalPoint &outer) const
Definition: TrackFitter.cc:162
T curvature(T InversePt, const edm::EventSetup &iSetup)
OutputIterator zip(InputIterator1 first1, InputIterator1 last1, InputIterator2 first2, InputIterator2 last2, OutputIterator result, Compare comp)
reco::Track * build(const Measurement1D &pt, const Measurement1D &phi, const Measurement1D &cotTheta, const Measurement1D &tip, const Measurement1D &zip, float chi2, int charge, const std::vector< const TrackingRecHit * > &hits, const MagneticField *mf, const GlobalPoint &reference=GlobalPoint(0, 0, 0)) const
float getCotThetaAndUpdateZip(const GlobalPoint &inner, const GlobalPoint &outer, float radius, float phi, float d0, float &zip) const
Definition: TrackFitter.cc:130
int getCharge(int phi1, int phi2, int phi3, int phi4, int mode)
T y() const
Cartesian y coordinate.
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
Definition: RZLine.h:12
const edm::EventSetup * theES
Definition: TrackFitter.h:38
Definition: errors.py:1
#define declareDynArray(T, n, x)
Definition: DynArray.h:59
float getPhi(float xC, float yC, int charge) const
Definition: TrackFitter.cc:150
const TransientTrackingRecHitBuilder * theTTRecHitBuilder
Definition: TrackFitter.h:41
T x() const
Cartesian x coordinate.
const MagneticField * theField
Definition: TrackFitter.h:40

Member Data Documentation

const edm::EventSetup* TrackFitter::theES
private

Definition at line 38 of file TrackFitter.h.

const MagneticField* TrackFitter::theField
private

Definition at line 40 of file TrackFitter.h.

const TrackerGeometry* TrackFitter::theTracker
private

Definition at line 39 of file TrackFitter.h.

const TransientTrackingRecHitBuilder* TrackFitter::theTTRecHitBuilder
private

Definition at line 41 of file TrackFitter.h.