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 TrackerGeometry *tracker, const MagneticField *field, const TransientTrackingRecHitBuilder *ttrhBuilder)
 
 ~TrackFitter () override
 
- Public Member Functions inherited from PixelFitterBase
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 MagneticFieldtheField
 
const TrackerGeometrytheTracker
 
const TransientTrackingRecHitBuildertheTTRecHitBuilder
 

Detailed Description

Definition at line 16 of file TrackFitter.h.

Constructor & Destructor Documentation

◆ TrackFitter()

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

Definition at line 18 of file TrackFitter.h.

21  : theTracker(tracker), theField(field), theTTRecHitBuilder(ttrhBuilder) {}
const TrackerGeometry * theTracker
Definition: TrackFitter.h:34
const TransientTrackingRecHitBuilder * theTTRecHitBuilder
Definition: TrackFitter.h:36
const MagneticField * theField
Definition: TrackFitter.h:35

◆ ~TrackFitter()

TrackFitter::~TrackFitter ( )
inlineoverride

Definition at line 22 of file TrackFitter.h.

22 {}

Member Function Documentation

◆ getCotThetaAndUpdateZip()

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

Definition at line 110 of file TrackFitter.cc.

References funct::cos(), d0, PVValHelper::dz, SurfaceOrientation::inner, listHistos::IP, M_PI_2, SurfaceOrientation::outer, perp(), CosmicsPD_Skims::radius, funct::sin(), and reco::zip().

111  {
112  float chi = phi - M_PI_2;
113  GlobalPoint IP(d0 * cos(chi), d0 * sin(chi), zip);
114 
115  float phi1 = 2 * asin(0.5 * (inner - IP).perp() / radius);
116  float phi2 = 2 * asin(0.5 * (outer - IP).perp() / radius);
117 
118  float dr = radius * (phi2 - phi1);
119  float dz = outer.z() - inner.z();
120 
121  // Recalculate ZIP
122  zip = (inner.z() * phi2 - outer.z() * phi1) / (phi2 - phi1);
123 
124  return (fabs(dr) > 1.e-3) ? dz / dr : 0;
125 }
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float zip(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:90
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define M_PI_2
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
T perp() const
Magnitude of transverse component.
static constexpr float d0

◆ getErrTipAndErrZip()

void TrackFitter::getErrTipAndErrZip ( float  pt,
float  eta,
float &  errZip,
float &  errTip 
) const
private

Definition at line 157 of file TrackFitter.cc.

References PVValHelper::eta, DiDispStaMuonMonitor_cfi::pt, and mathSSE::sqrt().

157  {
158  float coshEta = cosh(eta);
159 
160  { // transverse
161  float c_ms = 0.0115; //0.0115;
162  float s_le = 0.0095; //0.0123;
163  float s_ms2 = c_ms * c_ms / (pt * pt) * coshEta;
164 
165  errTip = sqrt(s_le * s_le + s_ms2);
166  }
167 
168  { // z
169  float c_ms = 0.0070;
170  float s_le = 0.0135;
171 
172  errZip = sqrt((s_le * s_le + c_ms * c_ms / (pt * pt)) * coshEta * coshEta * coshEta);
173  }
174 }
T sqrt(T t)
Definition: SSEVec.h:23

◆ getPhi()

float TrackFitter::getPhi ( float  xC,
float  yC,
int  charge 
) const
private

Definition at line 128 of file TrackFitter.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge.

128  {
129  float phiC;
130 
131  if (charge > 0)
132  phiC = atan2(xC, -yC);
133  else
134  phiC = atan2(-xC, yC);
135 
136  return phiC;
137 }

◆ getZip()

float TrackFitter::getZip ( float  d0,
float  curv,
const GlobalPoint inner,
const GlobalPoint outer 
) const
private

Definition at line 140 of file TrackFitter.cc.

References SurfaceOrientation::inner, SurfaceOrientation::outer, diffTwoXMLs::r2, ppsModifySingularModes_cfi::z1, and ppsModifySingularModes_cfi::z2.

140  {
141  // phi = asin(r*rho/2) with asin(x) ~= x+x**3/(2*3)
142  float rho3 = curv * curv * curv;
143 
144  float r1 = inner.perp();
145  double phi1 = r1 * curv / 2 + inner.perp2() * r1 * rho3 / 48.;
146 
147  float r2 = outer.perp();
148  double phi2 = r2 * curv / 2 + outer.perp2() * r2 * rho3 / 48.;
149 
150  double z1 = inner.z();
151  double z2 = outer.z();
152 
153  return z1 - phi1 / (phi1 - phi2) * (z1 - z2);
154 }

◆ run()

std::unique_ptr< reco::Track > TrackFitter::run ( const std::vector< const TrackingRecHit *> &  hits,
const TrackingRegion region 
) const
overridevirtual

Implements PixelFitterBase.

Definition at line 41 of file TrackFitter.cc.

References PixelTrackBuilder::build(), CircleFromThreePoints::center(), ALCARECOTkAlJpsiMuMu_cff::charge, isoTrack_cff::chi2, RZLine::chi2(), PixelRecoUtilities::curvature(), CircleFromThreePoints::curvature(), declareDynArray, vertexPlots::e4, deDxTools::getCharge(), mkfit::getPhi(), hfClusterShapes_cfi::hits, mps_fire::i, PixelRecoUtilities::inversePt(), PixelPluginsPhase0_cfi::isBarrel, Basic2DVector< T >::mag(), TrackingDataMCValidation_Standalone_cff::nhits, DiDispStaMuonMonitor_cfi::pt, rpcPointValidation_cfi::recHit, runTheMatrix::ret, qcdUeDQM_cfi::tip, Basic2DVector< T >::x(), Basic2DVector< T >::y(), and reco::zip().

42  {
43  std::unique_ptr<reco::Track> ret;
44 
45  int nhits = hits.size();
46  if (nhits < 2)
47  return ret;
48 
52 
53  unsigned int i = 0;
54  for (auto const& ih : hits) {
55  auto recHit = theTTRecHitBuilder->build(ih);
56 
57  points[i] = recHit->globalPosition();
58  errors[i] = recHit->globalPositionError();
59  isBarrel[++i] = recHit->detUnit()->type().isBarrel();
60  }
61 
62  CircleFromThreePoints circle = (nhits == 2) ? CircleFromThreePoints(GlobalPoint(0., 0., 0.), points[0], points[1])
63  : CircleFromThreePoints(points[0], points[1], points[2]);
64 
65  int charge = getCharge(points);
66  float curvature = circle.curvature();
67 
68  // pt
70  float valPt = (invPt > 1.e-4) ? 1. / invPt : 1.e4;
71  float errPt = 0.055 * valPt + 0.017 * valPt * valPt;
72 
73  CircleFromThreePoints::Vector2D center = circle.center();
74 
75  // tip
76  float valTip = charge * (center.mag() - 1 / curvature);
77  // zip
78  float valZip = getZip(valTip, curvature, points[0], points[1]);
79  // phi
80  float valPhi = getPhi(center.x(), center.y(), charge);
81  // cot(theta), update zip
82  float valCotTheta = getCotThetaAndUpdateZip(points[0], points[1], 1 / curvature, valPhi, valTip, valZip);
83 
84  // errors
85  float errTip, errZip;
86  getErrTipAndErrZip(valPt, points.back().eta(), errTip, errZip);
87  float errPhi = 0.002;
88  float errCotTheta = 0.002;
89 
90  float chi2 = 0;
91  if (nhits > 2) {
92  RZLine rzLine(points, errors, isBarrel);
93  chi2 = rzLine.chi2();
94  }
95 
96  // build pixel track
97  PixelTrackBuilder builder;
98 
99  Measurement1D pt(valPt, errPt);
100  Measurement1D phi(valPhi, errPhi);
101  Measurement1D cotTheta(valCotTheta, errCotTheta);
102  Measurement1D tip(valTip, errTip);
103  Measurement1D zip(valZip, errZip);
104 
105  ret.reset(builder.build(pt, phi, cotTheta, tip, zip, chi2, charge, hits, theField));
106  return ret;
107 }
ret
prodAgent to be discontinued
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float zip(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:90
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
float getZip(float d0, float curv, const GlobalPoint &inner, const GlobalPoint &outer) const
Definition: TrackFitter.cc:140
T y() const
Cartesian y coordinate.
int getCharge(const SiStripCluster *cluster, int &nSatStrip, const GeomDetUnit &detUnit, const std::vector< std::vector< float >> &calibGains, const unsigned int &m_off)
float getPhi(float xC, float yC, int charge) const
Definition: TrackFitter.cc:128
T curvature(T InversePt, const MagneticField &field)
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
void getErrTipAndErrZip(float pt, float eta, float &errZip, float &errTip) const
Definition: TrackFitter.cc:157
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
T inversePt(T curvature, const MagneticField &field)
Definition: RZLine.h:12
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
float getCotThetaAndUpdateZip(const GlobalPoint &inner, const GlobalPoint &outer, float radius, float phi, float d0, float &zip) const
Definition: TrackFitter.cc:110
Definition: errors.py:1
#define declareDynArray(T, n, x)
Definition: DynArray.h:91
const TransientTrackingRecHitBuilder * theTTRecHitBuilder
Definition: TrackFitter.h:36
T x() const
Cartesian x coordinate.
const MagneticField * theField
Definition: TrackFitter.h:35

Member Data Documentation

◆ theField

const MagneticField* TrackFitter::theField
private

Definition at line 35 of file TrackFitter.h.

◆ theTracker

const TrackerGeometry* TrackFitter::theTracker
private

Definition at line 34 of file TrackFitter.h.

◆ theTTRecHitBuilder

const TransientTrackingRecHitBuilder* TrackFitter::theTTRecHitBuilder
private

Definition at line 36 of file TrackFitter.h.