PandaRoot
PndFTSCATrackParamVector Class Reference

#include <PndFTSCATrackParamVector.h>

Classes

struct  PndFTSCATrackFitParam
 

Public Member Functions

 PndFTSCATrackParamVector ()
 
void ConvertTrackParamToVector (PndFTSCATrackParam t0[float_v::Size], int nTracksV)
 
void InitCovMatrix (float_v d2QMom=0.f)
 
void InitByTarget (const FTSCATarget &target)
 
void InitByHit (const FTSCAHitV &hit, const PndFTSCAParam &param, const float_v &dQMom)
 
void InitDirection (float_v r0, float_v r1, float_v r2)
 
float_v X () const
 
float_v Y () const
 
float_v Z () const
 
float_v SinPhi () const
 
float_v DzDs () const
 
float_v QPt () const
 
float_v Angle () const
 
float_v X0 () const
 
float_v X1 () const
 
float_v X2 () const
 
float_v Tx1 () const
 
float_v Tx2 () const
 
float_v QMomentum () const
 
float_v SignCosPhi () const
 
float_v Chi2 () const
 
int_v NDF () const
 
float_v Err2Y () const
 
float_v Err2Z () const
 
float_v Err2SinPhi () const
 
float_v Err2DzDs () const
 
float_v Err2QPt () const
 
float_v GetX () const
 
float_v GetY () const
 
float_v GetZ () const
 
float_v GetSinPhi () const
 
float_v GetDzDs () const
 
float_v GetQPt () const
 
float_v GetSignCosPhi () const
 
float_v GetChi2 () const
 
int_v GetNDF () const
 
float_v GetKappa (const float_v &Bz) const
 
float_v GetCosPhiPositive () const
 
float_v GetCosPhi () const
 
float_v Err2X1 () const
 
float_v Err2X2 () const
 
float_v Err2QMomentum () const
 
const float_v & Par (int i) const
 
const float_v & Cov (int i) const
 
void SetTrackParam (const PndFTSCATrackParamVector &param, const float_m &m=float_m(true))
 
void SetTrackParamOne (int iV, const PndFTSCATrackParamVector &param, int iVa)
 
void SetPar (int i, const float_v &v)
 
void SetPar (int i, const float_v &v, const float_m &m)
 
void SetCov (int i, const float_v &v)
 
void SetCov (int i, const float_v &v, const float_m &m)
 
void SetX (const float_v &v)
 
void SetY (const float_v &v)
 
void SetZ (const float_v &v)
 
void SetX (const float_v &v, const float_m &m)
 
void SetY (const float_v &v, const float_m &m)
 
void SetZ (const float_v &v, const float_m &m)
 
void SetSinPhi (const float_v &v)
 
void SetSinPhi (const float_v &v, const float_m &m)
 
void SetDzDs (const float_v &v)
 
void SetDzDs (const float_v &v, const float_m &m)
 
void SetQPt (const float_v &v)
 
void SetQMomentum (const float_v &v)
 
void SetQPt (const float_v &v, const float_m &m)
 
void SetSignCosPhi (const float_v &v)
 
void SetSignCosPhi (const float_v &v, const float_m &m)
 
void SetChi2 (const float_v &v)
 
void SetChi2 (const float_v &v, const float_m &m)
 
void SetNDF (int v)
 
void SetNDF (const int_v &v)
 
void SetNDF (const int_v &v, const int_m &m)
 
void SetAngle (const float_v &v)
 
void SetAngle (const float_v &v, const float_m &m)
 
void SetErr2Y (float_v v)
 
void SetErr2Z (float_v v)
 
void SetErr2QPt (float_v v)
 
float_v GetDist2 (const PndFTSCATrackParamVector &t) const
 
float_v GetDistXZ2 (const PndFTSCATrackParamVector &t) const
 
float_v GetS (const float_v &x, const float_v &y, const float_v &Bz) const
 
void GetDCAPoint (const float_v &x, const float_v &y, const float_v &z, float_v *px, float_v *py, float_v *pz, const float_v &Bz) const
 
float_m TransportToX0WithMaterial (const float_v &x, const float_v &XOverX0, const float_v &XThimesRho, const float_v &Bz, const float maxSinPhi=.999f)
 
float_m TransportToX0 (const float_v &x, const float_v &Bz, const float maxSinPhi=.999f, const float_m &mask=float_m(true))
 
float_m TransportToX0 (const float_v &x, PndFTSCATrackLinearisationVector &t0, const float_v &Bz, const float maxSinPhi=.999f, float_v *DL=0, const float_m &mask=float_m(true))
 
float_m TransportToX0 (const float_v &x, const float_v &sinPhi0, const float_v &Bz, const float_v maxSinPhi=.999f, const float_m &mask=float_m(true))
 
float_m TransportToX0WithMaterial (const float_v &x, PndFTSCATrackLinearisationVector &t0, PndFTSCATrackFitParam &par, const float_v &XOverX0, const float_v &XThimesRho, const float_v &Bz, const float maxSinPhi=.999f, const float_m &mask=float_m(true))
 
float_m TransportToX0WithMaterial (const float_v &x, PndFTSCATrackFitParam &par, const float_v &XOverX0, const float_v &XThimesRho, const float_v &Bz, const float maxSinPhi=.999f)
 
float_m Rotate (const float_v &alpha, PndFTSCATrackLinearisationVector &t0, const float maxSinPhi=.999f, const float_m &mask=float_m(true))
 
float_m Rotate (const float_v &alpha, const float maxSinPhi=.999f, const float_m &mask=float_m(true))
 
void RotateXY (float_v alpha, float_v &x, float_v &y, float_v &sin, const float_m &mask=float_m(true)) const
 
float_m FilterWithMaterial (const float_v &y, const float_v &z, float_v err2Y, float_v errYZ, float_v err2Z, float maxSinPhi=0.999f, const float_m &mask=float_m(true), const int_v &hitNDF=int_v(2), const float_v &chi2Cut=10e10f)
 
float_m FilterWithMaterial (const float_v &y, const float_v &z, const FTSCAStripInfo &info, float_v err2, float maxSinPhi=0.999f, const float_m &mask=float_m(true), const float_v &chi2Cut=10e10f)
 
float_m FilterWithMaterial (const float_v &y, const float_v &z, const float_v &r, const FTSCAStripInfo &info, float_v err2, float maxSinPhi=0.999f, const float_m &mask=float_m(true), const float_v &chi2Cut=10e10f)
 
void CalculateFitParameters (PndFTSCATrackFitParam &par, const float_v &mass=0.13957f)
 
float_m CorrectForMeanMaterial (const float_v &xOverX0, const float_v &xTimesRho, const PndFTSCATrackFitParam &par, const float_m &_mask)
 
float_m FilterVtx (const float_v &xV, const float_v &yV, const CAX1X2MeasurementInfo &info, float_v &extrDx, float_v &extrDy, float_v J[], const float_m &active=float_m(true))
 
float_m TransportJ0ToX0 (const float_v &x0, const float_v &cBz, float_v &extrDy, float_v &extrDz, float_v J[], const float_m &active=float_m(true))
 
float_m Transport (const int_v &ista, const PndFTSCAParam &param, const float_m &mask=float_m(true))
 
float_m Transport (const FTSCAHitV &hit, const PndFTSCAParam &p, const float_m &mask=float_m(true))
 
float_m Filter (const FTSCAHitV &hit, const PndFTSCAParam &param, const float_m &mask=float_m(true), const float_v &chi2Cut=10e10f)
 
float_m Transport (const FTSCAHit &hit, const PndFTSCAParam &p, const float_m &mask=float_m(true))
 
float_m Filter (const FTSCAHit &hit, const PndFTSCAParam &param, const float_m &mask=float_m(true), const float_v &chi2Cut=10e10f)
 
float_m AddTarget (const FTSCATarget &target, const float_m &mask=float_m(true))
 

Static Public Member Functions

static float_v ApproximateBetheBloch (const float_v &beta2)
 
static float_v BetheBlochGeant (const float_v &bg, const float_v &kp0=2.33f, const float_v &kp1=0.20f, const float_v &kp2=3.00f, const float_v &kp3=173e-9f, const float_v &kp4=0.49848f)
 
static float_v BetheBlochSolid (const float_v &bg)
 
static float_v BetheBlochGas (const float_v &bg)
 

Friends

class PndFTSCATrackParam
 
std::istream & operator>> (std::istream &, PndFTSCATrackParamVector &)
 
std::ostream & operator<< (std::ostream &, const PndFTSCATrackParamVector &)
 

Detailed Description

PndFTSCATrackParamVector class describes the track parametrisation which is used by the PndFTSCATracker slice tracker.

Definition at line 1285 of file PndFTSCATrackParamVector.h.

Constructor & Destructor Documentation

◆ PndFTSCATrackParamVector()

PndFTSCATrackParamVector::PndFTSCATrackParamVector ( )
inline

Definition at line 1290 of file PndFTSCATrackParamVector.h.

References ConvertTrackParamToVector(), i, InitByHit(), InitByTarget(), and InitCovMatrix().

1290  : fX(Vc::Zero), fSignCosPhi(Vc::Zero), fChi2(Vc::Zero), fNDF(Vc::Zero)
1291  {
1292  for (int i = 0; i < 5; ++i)
1293  fP[i].setZero();
1294  for (int i = 0; i < 15; ++i)
1295  fC[i].setZero();
1296  }
unsigned int i
Definition: P4_F32vec4.h:33

Member Function Documentation

◆ AddTarget()

float_m PndFTSCATrackParamVector::AddTarget ( const FTSCATarget target,
const float_m &  mask = float_m(true) 
)

Referenced by SetErr2QPt().

◆ Angle()

float_v PndFTSCATrackParamVector::Angle ( ) const
inline

Definition at line 1331 of file PndFTSCATrackParamVector.h.

Referenced by SetTrackParam(), and SetTrackParamOne().

1331 { return fAlpha; }

◆ ApproximateBetheBloch()

static float_v PndFTSCATrackParamVector::ApproximateBetheBloch ( const float_v &  beta2)
static

Referenced by SetErr2QPt().

◆ BetheBlochGas()

static float_v PndFTSCATrackParamVector::BetheBlochGas ( const float_v &  bg)
static

Referenced by SetErr2QPt().

◆ BetheBlochGeant()

static float_v PndFTSCATrackParamVector::BetheBlochGeant ( const float_v &  bg,
const float_v &  kp0 = 2.33f,
const float_v &  kp1 = 0.20f,
const float_v &  kp2 = 3.00f,
const float_v &  kp3 = 173e-9f,
const float_v &  kp4 = 0.49848f 
)
static

Referenced by SetErr2QPt().

◆ BetheBlochSolid()

static float_v PndFTSCATrackParamVector::BetheBlochSolid ( const float_v &  bg)
static

Referenced by SetErr2QPt().

◆ CalculateFitParameters()

void PndFTSCATrackParamVector::CalculateFitParameters ( PndFTSCATrackFitParam par,
const float_v &  mass = 0.13957f 
)

Referenced by SetErr2QPt().

◆ Chi2()

float_v PndFTSCATrackParamVector::Chi2 ( ) const
inline

Definition at line 1345 of file PndFTSCATrackParamVector.h.

Referenced by PndFTSCAGBTracker::Refit_1().

1345 { return fChi2; }

◆ ConvertTrackParamToVector()

void PndFTSCATrackParamVector::ConvertTrackParamToVector ( PndFTSCATrackParam  t0[float_v::Size],
int  nTracksV 
)

◆ CorrectForMeanMaterial()

float_m PndFTSCATrackParamVector::CorrectForMeanMaterial ( const float_v &  xOverX0,
const float_v &  xTimesRho,
const PndFTSCATrackFitParam par,
const float_m &  _mask 
)

Referenced by SetErr2QPt().

◆ Cov()

const float_v& PndFTSCATrackParamVector::Cov ( int  i) const
inline

Definition at line 1375 of file PndFTSCATrackParamVector.h.

References i.

Referenced by SetTrackParam(), and SetTrackParamOne().

1375 { return fC[i]; }
unsigned int i
Definition: P4_F32vec4.h:33

◆ DzDs()

float_v PndFTSCATrackParamVector::DzDs ( ) const
inline

Definition at line 1328 of file PndFTSCATrackParamVector.h.

Referenced by TransportJ0ToX0(), TransportToX0(), and Tx2().

1328 { return fP[3]; }

◆ Err2DzDs()

float_v PndFTSCATrackParamVector::Err2DzDs ( ) const
inline

Definition at line 1351 of file PndFTSCATrackParamVector.h.

1351 { return fC[9]; }

◆ Err2QMomentum()

float_v PndFTSCATrackParamVector::Err2QMomentum ( ) const
inline

Definition at line 1372 of file PndFTSCATrackParamVector.h.

1372 { return fC[14]; }

◆ Err2QPt()

float_v PndFTSCATrackParamVector::Err2QPt ( ) const
inline

Definition at line 1352 of file PndFTSCATrackParamVector.h.

1352 { return fC[14]; }

◆ Err2SinPhi()

float_v PndFTSCATrackParamVector::Err2SinPhi ( ) const
inline

Definition at line 1350 of file PndFTSCATrackParamVector.h.

1350 { return fC[5]; }

◆ Err2X1()

float_v PndFTSCATrackParamVector::Err2X1 ( ) const
inline

Definition at line 1368 of file PndFTSCATrackParamVector.h.

1368 { return fC[0]; }

◆ Err2X2()

float_v PndFTSCATrackParamVector::Err2X2 ( ) const
inline

Definition at line 1369 of file PndFTSCATrackParamVector.h.

1369 { return fC[2]; }

◆ Err2Y()

float_v PndFTSCATrackParamVector::Err2Y ( ) const
inline

Definition at line 1348 of file PndFTSCATrackParamVector.h.

1348 { return fC[0]; }

◆ Err2Z()

float_v PndFTSCATrackParamVector::Err2Z ( ) const
inline

Definition at line 1349 of file PndFTSCATrackParamVector.h.

1349 { return fC[2]; }

◆ Filter() [1/2]

float_m PndFTSCATrackParamVector::Filter ( const FTSCAHitV hit,
const PndFTSCAParam param,
const float_m &  mask = float_m(true),
const float_v &  chi2Cut = 10e10f 
)

Referenced by SetErr2QPt().

◆ Filter() [2/2]

float_m PndFTSCATrackParamVector::Filter ( const FTSCAHit hit,
const PndFTSCAParam param,
const float_m &  mask = float_m(true),
const float_v &  chi2Cut = 10e10f 
)

◆ FilterVtx()

float_m PndFTSCATrackParamVector::FilterVtx ( const float_v &  xV,
const float_v &  yV,
const CAX1X2MeasurementInfo info,
float_v &  extrDx,
float_v &  extrDy,
float_v  J[],
const float_m &  active = float_m(true) 
)
inline

Definition at line 1914 of file PndFTSCATrackParamVector.h.

References CAX1X2MeasurementInfo::C00(), CAX1X2MeasurementInfo::C10(), CAX1X2MeasurementInfo::C11(), Y(), and Z().

Referenced by SetErr2QPt().

1915 {
1916  float_v zeta0, zeta1, S00, S10, S11, si;
1917  float_v F00, F10, F20, F30, F40, F01, F11, F21, F31, F41;
1918  float_v K00, K10, K20, K30, K40, K01, K11, K21, K31, K41;
1919  float_v &c00 = fC[0];
1920  float_v &c10 = fC[1];
1921  float_v &c11 = fC[2];
1922  float_v &c20 = fC[3];
1923  float_v &c21 = fC[4];
1924  float_v &c22 = fC[5];
1925  float_v &c30 = fC[6];
1926  float_v &c31 = fC[7];
1927  float_v &c32 = fC[8];
1928  float_v &c33 = fC[9];
1929  float_v &c40 = fC[10];
1930  float_v &c41 = fC[11];
1931  float_v &c42 = fC[12];
1932  float_v &c43 = fC[13];
1933  float_v &c44 = fC[14];
1934 
1935  zeta0 = Y() + extrDy - yV;
1936  zeta1 = Z() + extrDz - zV;
1937 
1938  // H = 1 0 J[0] J[1] J[2]
1939  // 0 1 J[3] J[4] J[5]
1940 
1941  // F = CH'
1942  F00 = c00;
1943  F01 = c10;
1944  F10 = c10;
1945  F11 = c11;
1946  F20 = J[0] * c22;
1947  F21 = J[3] * c22;
1948  F30 = J[1] * c33;
1949  F31 = J[4] * c33;
1950  F40 = J[2] * c44;
1951  F41 = J[5] * c44;
1952 
1953  S00 = info.C00() + F00 + J[0] * F20 + J[1] * F30 + J[2] * F40;
1954  S10 = info.C10() + F10 + J[3] * F20 + J[4] * F30 + J[5] * F40;
1955  S11 = info.C11() + F11 + J[3] * F21 + J[4] * F31 + J[5] * F41;
1956 
1957  si = 1.f / (S00 * S11 - S10 * S10);
1958  float_v S00tmp = S00;
1959  S00 = si * S11;
1960  S10 = -si * S10;
1961  S11 = si * S00tmp;
1962 
1963  fChi2(active) += zeta0 * zeta0 * S00 + 2.f * zeta0 * zeta1 * S10 + zeta1 * zeta1 * S11;
1964 
1965  K00 = F00 * S00 + F01 * S10;
1966  K01 = F00 * S10 + F01 * S11;
1967  K10 = F10 * S00 + F11 * S10;
1968  K11 = F10 * S10 + F11 * S11;
1969  K20 = F20 * S00 + F21 * S10;
1970  K21 = F20 * S10 + F21 * S11;
1971  K30 = F30 * S00 + F31 * S10;
1972  K31 = F30 * S10 + F31 * S11;
1973  K40 = F40 * S00 + F41 * S10;
1974  K41 = F40 * S10 + F41 * S11;
1975 
1976  fP[0](active) -= K00 * zeta0 + K01 * zeta1;
1977  fP[1](active) -= K10 * zeta0 + K11 * zeta1;
1978  fP[2](active) -= K20 * zeta0 + K21 * zeta1;
1979  fP[3](active) -= K30 * zeta0 + K31 * zeta1;
1980  fP[4](active) -= K40 * zeta0 + K41 * zeta1;
1981 
1982  c00(active) -= (K00 * F00 + K01 * F01);
1983  c10(active) -= (K10 * F00 + K11 * F01);
1984  c11(active) -= (K10 * F10 + K11 * F11);
1985  c20(active) = -(K20 * F00 + K21 * F01);
1986  c21(active) = -(K20 * F10 + K21 * F11);
1987  c22(active) -= (K20 * F20 + K21 * F21);
1988  c30(active) = -(K30 * F00 + K31 * F01);
1989  c31(active) = -(K30 * F10 + K31 * F11);
1990  c32(active) = -(K30 * F20 + K31 * F21);
1991  c33(active) -= (K30 * F30 + K31 * F31);
1992  c40(active) = -(K40 * F00 + K41 * F01);
1993  c41(active) = -(K40 * F10 + K41 * F11);
1994  c42(active) = -(K40 * F20 + K41 * F21);
1995  c43(active) = -(K40 * F30 + K41 * F31);
1996  c44(active) -= (K40 * F40 + K41 * F41);
1997 
1998  return active;
1999 }
const float_v & C00() const
const float_v & C10() const
const float_v & C11() const

◆ FilterWithMaterial() [1/3]

float_m PndFTSCATrackParamVector::FilterWithMaterial ( const float_v &  y,
const float_v &  z,
float_v  err2Y,
float_v  errYZ,
float_v  err2Z,
float  maxSinPhi = 0.999f,
const float_m &  mask = float_m(true),
const int_v &  hitNDF = int_v(2),
const float_v &  chi2Cut = 10e10f 
)

Referenced by SetErr2QPt().

◆ FilterWithMaterial() [2/3]

float_m PndFTSCATrackParamVector::FilterWithMaterial ( const float_v &  y,
const float_v &  z,
const FTSCAStripInfo info,
float_v  err2,
float  maxSinPhi = 0.999f,
const float_m &  mask = float_m(true),
const float_v &  chi2Cut = 10e10f 
)

◆ FilterWithMaterial() [3/3]

float_m PndFTSCATrackParamVector::FilterWithMaterial ( const float_v &  y,
const float_v &  z,
const float_v &  r,
const FTSCAStripInfo info,
float_v  err2,
float  maxSinPhi = 0.999f,
const float_m &  mask = float_m(true),
const float_v &  chi2Cut = 10e10f 
)

◆ GetChi2()

float_v PndFTSCATrackParamVector::GetChi2 ( ) const
inline

Definition at line 1361 of file PndFTSCATrackParamVector.h.

Referenced by SetTrackParam(), and SetTrackParamOne().

1361 { return fChi2; }

◆ GetCosPhi()

float_v PndFTSCATrackParamVector::GetCosPhi ( ) const
inline

Definition at line 1366 of file PndFTSCATrackParamVector.h.

References SinPhi(), and CAMath::Sqrt().

Referenced by Rotate(), and RotateXY().

1366 { return fSignCosPhi * CAMath::Sqrt(float_v(Vc::One) - SinPhi() * SinPhi()); }
static T Sqrt(const T &x)
Definition: PndCAMath.h:57

◆ GetCosPhiPositive()

float_v PndFTSCATrackParamVector::GetCosPhiPositive ( ) const
inline

Definition at line 1365 of file PndFTSCATrackParamVector.h.

References SinPhi(), and CAMath::Sqrt().

1365 { return CAMath::Sqrt(float_v(Vc::One) - SinPhi() * SinPhi()); }
static T Sqrt(const T &x)
Definition: PndCAMath.h:57

◆ GetDCAPoint()

void PndFTSCATrackParamVector::GetDCAPoint ( const float_v &  x,
const float_v &  y,
const float_v &  z,
float_v *  px,
float_v *  py,
float_v *  pz,
const float_v &  Bz 
) const

Referenced by SetErr2QPt().

◆ GetDist2()

float_v PndFTSCATrackParamVector::GetDist2 ( const PndFTSCATrackParamVector t) const

Referenced by SetErr2QPt().

◆ GetDistXZ2()

float_v PndFTSCATrackParamVector::GetDistXZ2 ( const PndFTSCATrackParamVector t) const

Referenced by SetErr2QPt().

◆ GetDzDs()

float_v PndFTSCATrackParamVector::GetDzDs ( ) const
inline

Definition at line 1358 of file PndFTSCATrackParamVector.h.

1358 { return fP[3]; }

◆ GetKappa()

float_v PndFTSCATrackParamVector::GetKappa ( const float_v &  Bz) const
inline

Definition at line 1364 of file PndFTSCATrackParamVector.h.

1364 { return fP[4] * Bz; }

◆ GetNDF()

int_v PndFTSCATrackParamVector::GetNDF ( ) const
inline

Definition at line 1362 of file PndFTSCATrackParamVector.h.

Referenced by SetTrackParam(), and SetTrackParamOne().

1362 { return fNDF; }

◆ GetQPt()

float_v PndFTSCATrackParamVector::GetQPt ( ) const
inline

Definition at line 1359 of file PndFTSCATrackParamVector.h.

1359 { return fP[4]; }

◆ GetS()

float_v PndFTSCATrackParamVector::GetS ( const float_v &  x,
const float_v &  y,
const float_v &  Bz 
) const

Referenced by SetErr2QPt().

◆ GetSignCosPhi()

float_v PndFTSCATrackParamVector::GetSignCosPhi ( ) const
inline

Definition at line 1360 of file PndFTSCATrackParamVector.h.

1360 { return fSignCosPhi; }

◆ GetSinPhi()

float_v PndFTSCATrackParamVector::GetSinPhi ( ) const
inline

Definition at line 1357 of file PndFTSCATrackParamVector.h.

1357 { return fP[2]; }

◆ GetX()

float_v PndFTSCATrackParamVector::GetX ( ) const
inline

Definition at line 1354 of file PndFTSCATrackParamVector.h.

1354 { return fX; }

◆ GetY()

float_v PndFTSCATrackParamVector::GetY ( ) const
inline

Definition at line 1355 of file PndFTSCATrackParamVector.h.

1355 { return fP[0]; }

◆ GetZ()

float_v PndFTSCATrackParamVector::GetZ ( ) const
inline

Definition at line 1356 of file PndFTSCATrackParamVector.h.

1356 { return fP[1]; }

◆ InitByHit()

void PndFTSCATrackParamVector::InitByHit ( const FTSCAHitV hit,
const PndFTSCAParam param,
const float_v &  dQMom 
)

◆ InitByTarget()

void PndFTSCATrackParamVector::InitByTarget ( const FTSCATarget target)

◆ InitCovMatrix()

void PndFTSCATrackParamVector::InitCovMatrix ( float_v  d2QMom = 0.f)

◆ InitDirection()

void PndFTSCATrackParamVector::InitDirection ( float_v  r0,
float_v  r1,
float_v  r2 
)
inline

Definition at line 1304 of file PndFTSCATrackParamVector.h.

References SetDzDs(), SetSignCosPhi(), SetSinPhi(), and sqrt().

1305  {
1306  const float_v r = sqrt(r0 * r0 + r1 * r1);
1307  SetSinPhi(r1 / r);
1308  SetSignCosPhi(r0 / abs(r0));
1309  SetDzDs(r2 / r);
1310  }
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:40
void SetSignCosPhi(const float_v &v)
void SetSinPhi(const float_v &v)

◆ NDF()

int_v PndFTSCATrackParamVector::NDF ( ) const
inline

Definition at line 1346 of file PndFTSCATrackParamVector.h.

Referenced by PndFTSCAGBTracker::Refit_1().

1346 { return fNDF; }

◆ Par()

const float_v& PndFTSCATrackParamVector::Par ( int  i) const
inline

Definition at line 1374 of file PndFTSCATrackParamVector.h.

References i.

Referenced by SetTrackParam(), and SetTrackParamOne().

1374 { return fP[i]; }
unsigned int i
Definition: P4_F32vec4.h:33

◆ QMomentum()

float_v PndFTSCATrackParamVector::QMomentum ( ) const
inline

Definition at line 1338 of file PndFTSCATrackParamVector.h.

References QPt().

1338 { return QPt(); } // used for triplets comparison

◆ QPt()

float_v PndFTSCATrackParamVector::QPt ( ) const
inline

Definition at line 1329 of file PndFTSCATrackParamVector.h.

Referenced by QMomentum(), and TransportToX0().

1329 { return fP[4]; }

◆ Rotate() [1/2]

float_m PndFTSCATrackParamVector::Rotate ( const float_v &  alpha,
PndFTSCATrackLinearisationVector t0,
const float  maxSinPhi = .999f,
const float_m &  mask = float_m(true) 
)

Referenced by SetErr2QPt().

◆ Rotate() [2/2]

float_m PndFTSCATrackParamVector::Rotate ( const float_v &  alpha,
const float  maxSinPhi = .999f,
const float_m &  mask = float_m(true) 
)
inline

Definition at line 1851 of file PndFTSCATrackParamVector.h.

References CAMath::Abs(), alpha, CAMath::Cos(), f, GetCosPhi(), SetSignCosPhi(), SetSinPhi(), SetX(), SetY(), CAMath::Sin(), SinPhi(), X(), and Y().

1852 {
1853  //* Rotate the coordinate system in XY on the angle alpha
1854  if ((abs(alpha) < 1e-6f || !mask).isFull())
1855  return mask;
1856 
1857  const float_v cA = CAMath::Cos(alpha);
1858  const float_v sA = CAMath::Sin(alpha);
1859  const float_v x = X(), y = Y(), sP = SinPhi(), cP = GetCosPhi();
1860  const float_v cosPhi = cP * cA + sP * sA;
1861  const float_v sinPhi = -cP * sA + sP * cA;
1862 
1863  float_m mReturn = mask && (CAMath::Abs(sinPhi) < maxSinPhi) && (CAMath::Abs(cosPhi) > 1.e-2f) && (CAMath::Abs(cP) > 1.e-2f);
1864 
1865  mReturn &= abs(alpha) < 3.1415f * 0.25f; // allow turn by 45 degree only
1866 
1867  const float_v j0 = cP / cosPhi;
1868  const float_v j2 = cosPhi / cP;
1869 
1870  SetX(x * cA + y * sA, mReturn);
1871  SetY(-x * sA + y * cA, mReturn);
1872  SetSignCosPhi(CAMath::Abs(cosPhi) / cosPhi, mReturn);
1873  SetSinPhi(sinPhi, mReturn);
1874 
1875  // float J[5][5] = { { j0, 0, 0, 0, 0 }, // Y
1876  // { 0, 1, 0, 0, 0 }, // Z
1877  // { 0, 0, j2, 0, 0 }, // SinPhi
1878  // { 0, 0, 0, 1, 0 }, // DzDs
1879  // { 0, 0, 0, 0, 1 } }; // Kappa
1880  // cout<<"alpha="<<alpha<<" "<<x<<" "<<y<<" "<<sP<<" "<<cP<<" "<<j0<<" "<<j2<<endl;
1881  // cout<<" "<<fC[0]<<" "<<fC[1]<<" "<<fC[6]<<" "<<fC[10]<<" "<<fC[4]<<" "<<fC[5]<<" "<<fC[8]<<" "<<fC[12]<<endl;
1882  fC[0](mReturn) *= j0 * j0;
1883  fC[1](mReturn) *= j0;
1884  fC[3](mReturn) *= j0;
1885  fC[6](mReturn) *= j0;
1886  fC[10](mReturn) *= j0;
1887 
1888  fC[3](mReturn) *= j2;
1889  fC[4](mReturn) *= j2;
1890  fC[5](mReturn) *= j2 * j2;
1891  fC[8](mReturn) *= j2;
1892  fC[12](mReturn) *= j2;
1893 
1894  fAlpha(mReturn) += alpha;
1895  // cout<<" "<<fC[0]<<" "<<fC[1]<<" "<<fC[6]<<" "<<fC[10]<<" "<<fC[4]<<" "<<fC[5]<<" "<<fC[8]<<" "<<fC[12]<<endl;
1896  return mReturn;
1897 }
void SetSignCosPhi(const float_v &v)
static T Sin(const T &x)
Definition: PndCAMath.h:83
void SetSinPhi(const float_v &v)
static T Cos(const T &x)
Definition: PndCAMath.h:88
static T Abs(const T &x)
Definition: PndCAMath.h:68
float f
Definition: P4_F32vec4.h:32
double alpha
Definition: f_Init.h:31

◆ RotateXY()

void PndFTSCATrackParamVector::RotateXY ( float_v  alpha,
float_v &  x,
float_v &  y,
float_v &  sin,
const float_m &  mask = float_m(true) 
) const
inline

Definition at line 1899 of file PndFTSCATrackParamVector.h.

References CAMath::Cos(), GetCosPhi(), CAMath::Sin(), sin(), SinPhi(), X(), and Y().

Referenced by SetErr2QPt().

1900 {
1901  //* Rotate the coordinate system in XY on the angle alpha
1902  if ((abs(alpha) < 1e-6f || !mask).isFull())
1903  return;
1904 
1905  const float_v cA = CAMath::Cos(alpha);
1906  const float_v sA = CAMath::Sin(alpha);
1907 
1908  x(mask) = (X() * cA + Y() * sA);
1909  y(mask) = (-X() * sA + Y() * cA);
1910  sin(mask) = -GetCosPhi() * sA + SinPhi() * cA;
1911 }
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:130
static T Sin(const T &x)
Definition: PndCAMath.h:83
static T Cos(const T &x)
Definition: PndCAMath.h:88
float f
Definition: P4_F32vec4.h:32
double alpha
Definition: f_Init.h:31

◆ SetAngle() [1/2]

void PndFTSCATrackParamVector::SetAngle ( const float_v &  v)
inline

Definition at line 1437 of file PndFTSCATrackParamVector.h.

References v.

1437 { fAlpha = v; }
__m128 v
Definition: P4_F32vec4.h:15

◆ SetAngle() [2/2]

void PndFTSCATrackParamVector::SetAngle ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1438 of file PndFTSCATrackParamVector.h.

References v.

1438 { fAlpha(m) = v; }
__m128 m
Definition: P4_F32vec4.h:38
__m128 v
Definition: P4_F32vec4.h:15

◆ SetChi2() [1/2]

void PndFTSCATrackParamVector::SetChi2 ( const float_v &  v)
inline

Definition at line 1431 of file PndFTSCATrackParamVector.h.

References v.

Referenced by FTSCATrack::Fit2Times(), and PndFTSCAGBTracker::Refit_1().

1431 { fChi2 = v; }
__m128 v
Definition: P4_F32vec4.h:15

◆ SetChi2() [2/2]

void PndFTSCATrackParamVector::SetChi2 ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1432 of file PndFTSCATrackParamVector.h.

References v.

1432 { fChi2(m) = v; }
__m128 m
Definition: P4_F32vec4.h:38
__m128 v
Definition: P4_F32vec4.h:15

◆ SetCov() [1/2]

void PndFTSCATrackParamVector::SetCov ( int  i,
const float_v &  v 
)
inline

Definition at line 1413 of file PndFTSCATrackParamVector.h.

References i, and v.

1413 { fC[i] = v; }
__m128 v
Definition: P4_F32vec4.h:15
unsigned int i
Definition: P4_F32vec4.h:33

◆ SetCov() [2/2]

void PndFTSCATrackParamVector::SetCov ( int  i,
const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1414 of file PndFTSCATrackParamVector.h.

References i, and m.

1414 { fC[i](m) = v; }
__m128 m
Definition: P4_F32vec4.h:38
__m128 v
Definition: P4_F32vec4.h:15
unsigned int i
Definition: P4_F32vec4.h:33

◆ SetDzDs() [1/2]

void PndFTSCATrackParamVector::SetDzDs ( const float_v &  v)
inline

Definition at line 1424 of file PndFTSCATrackParamVector.h.

References v.

Referenced by InitDirection().

1424 { fP[3] = v; }
__m128 v
Definition: P4_F32vec4.h:15

◆ SetDzDs() [2/2]

void PndFTSCATrackParamVector::SetDzDs ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1425 of file PndFTSCATrackParamVector.h.

References m.

1425 { fP[3](m) = v; }
__m128 m
Definition: P4_F32vec4.h:38
__m128 v
Definition: P4_F32vec4.h:15

◆ SetErr2QPt()

◆ SetErr2Y()

void PndFTSCATrackParamVector::SetErr2Y ( float_v  v)
inline

Definition at line 1440 of file PndFTSCATrackParamVector.h.

References v.

1440 { fC[0] = v; }
__m128 v
Definition: P4_F32vec4.h:15

◆ SetErr2Z()

void PndFTSCATrackParamVector::SetErr2Z ( float_v  v)
inline

Definition at line 1441 of file PndFTSCATrackParamVector.h.

References v.

1441 { fC[2] = v; }
__m128 v
Definition: P4_F32vec4.h:15

◆ SetNDF() [1/3]

void PndFTSCATrackParamVector::SetNDF ( int  v)
inline

Definition at line 1433 of file PndFTSCATrackParamVector.h.

References v.

Referenced by FTSCATrack::Fit2Times().

1433 { fNDF = v; }
__m128 v
Definition: P4_F32vec4.h:15

◆ SetNDF() [2/3]

void PndFTSCATrackParamVector::SetNDF ( const int_v &  v)
inline

Definition at line 1434 of file PndFTSCATrackParamVector.h.

References v.

1434 { fNDF = v; }
__m128 v
Definition: P4_F32vec4.h:15

◆ SetNDF() [3/3]

void PndFTSCATrackParamVector::SetNDF ( const int_v &  v,
const int_m &  m 
)
inline

Definition at line 1435 of file PndFTSCATrackParamVector.h.

References v.

1435 { fNDF(m) = v; }
__m128 m
Definition: P4_F32vec4.h:38
__m128 v
Definition: P4_F32vec4.h:15

◆ SetPar() [1/2]

void PndFTSCATrackParamVector::SetPar ( int  i,
const float_v &  v 
)
inline

Definition at line 1411 of file PndFTSCATrackParamVector.h.

References i, and v.

1411 { fP[i] = v; }
__m128 v
Definition: P4_F32vec4.h:15
unsigned int i
Definition: P4_F32vec4.h:33

◆ SetPar() [2/2]

void PndFTSCATrackParamVector::SetPar ( int  i,
const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1412 of file PndFTSCATrackParamVector.h.

References i, and m.

1412 { fP[i](m) = v; }
__m128 m
Definition: P4_F32vec4.h:38
__m128 v
Definition: P4_F32vec4.h:15
unsigned int i
Definition: P4_F32vec4.h:33

◆ SetQMomentum()

void PndFTSCATrackParamVector::SetQMomentum ( const float_v &  v)
inline

Definition at line 1427 of file PndFTSCATrackParamVector.h.

References SetQPt().

1427 { SetQPt(v); }
__m128 v
Definition: P4_F32vec4.h:15

◆ SetQPt() [1/2]

void PndFTSCATrackParamVector::SetQPt ( const float_v &  v)
inline

Definition at line 1426 of file PndFTSCATrackParamVector.h.

References v.

Referenced by SetQMomentum().

1426 { fP[4] = v; }
__m128 v
Definition: P4_F32vec4.h:15

◆ SetQPt() [2/2]

void PndFTSCATrackParamVector::SetQPt ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1428 of file PndFTSCATrackParamVector.h.

References m.

1428 { fP[4](m) = v; }
__m128 m
Definition: P4_F32vec4.h:38
__m128 v
Definition: P4_F32vec4.h:15

◆ SetSignCosPhi() [1/2]

void PndFTSCATrackParamVector::SetSignCosPhi ( const float_v &  v)
inline

Definition at line 1429 of file PndFTSCATrackParamVector.h.

References v.

Referenced by InitDirection(), and Rotate().

1429 { fSignCosPhi = v; }
__m128 v
Definition: P4_F32vec4.h:15

◆ SetSignCosPhi() [2/2]

void PndFTSCATrackParamVector::SetSignCosPhi ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1430 of file PndFTSCATrackParamVector.h.

References v.

1430 { fSignCosPhi(m) = v; }
__m128 m
Definition: P4_F32vec4.h:38
__m128 v
Definition: P4_F32vec4.h:15

◆ SetSinPhi() [1/2]

void PndFTSCATrackParamVector::SetSinPhi ( const float_v &  v)
inline

Definition at line 1422 of file PndFTSCATrackParamVector.h.

References v.

Referenced by InitDirection(), and Rotate().

1422 { fP[2] = v; }
__m128 v
Definition: P4_F32vec4.h:15

◆ SetSinPhi() [2/2]

void PndFTSCATrackParamVector::SetSinPhi ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1423 of file PndFTSCATrackParamVector.h.

References m.

1423 { fP[2](m) = v; }
__m128 m
Definition: P4_F32vec4.h:38
__m128 v
Definition: P4_F32vec4.h:15

◆ SetTrackParam()

void PndFTSCATrackParamVector::SetTrackParam ( const PndFTSCATrackParamVector param,
const float_m &  m = float_m(true) 
)
inline

Definition at line 1385 of file PndFTSCATrackParamVector.h.

References Angle(), Cov(), GetChi2(), GetNDF(), i, m, Par(), SignCosPhi(), and X().

1386  {
1387  for (int i = 0; i < 5; i++)
1388  fP[i](m) = param.Par()[i];
1389  for (int i = 0; i < 15; i++)
1390  fC[i](m) = param.Cov()[i];
1391  fX(m) = param.X();
1392  fSignCosPhi(m) = param.SignCosPhi();
1393  fChi2(m) = param.GetChi2();
1394  fNDF(static_cast<int_m>(m)) = param.GetNDF();
1395  fAlpha(static_cast<int_m>(m)) = param.Angle();
1396  }
__m128 m
Definition: P4_F32vec4.h:38
unsigned int i
Definition: P4_F32vec4.h:33
const float_v & Par(int i) const
const float_v & Cov(int i) const

◆ SetTrackParamOne()

void PndFTSCATrackParamVector::SetTrackParamOne ( int  iV,
const PndFTSCATrackParamVector param,
int  iVa 
)
inline

Definition at line 1398 of file PndFTSCATrackParamVector.h.

References Angle(), Cov(), GetChi2(), GetNDF(), i, Par(), SignCosPhi(), and X().

Referenced by FTSCANPletV::CopyOne().

1399  {
1400  for (int i = 0; i < 5; i++)
1401  fP[i][iV] = param.Par()[i][iVa];
1402  for (int i = 0; i < 15; i++)
1403  fC[i][iV] = param.Cov()[i][iVa];
1404  fX[iV] = param.X()[iVa];
1405  fSignCosPhi[iV] = param.SignCosPhi()[iVa];
1406  fChi2[iV] = param.GetChi2()[iVa];
1407  fNDF[iV] = param.GetNDF()[iVa];
1408  fAlpha[iV] = param.Angle()[iVa];
1409  }
unsigned int i
Definition: P4_F32vec4.h:33
const float_v & Par(int i) const
const float_v & Cov(int i) const

◆ SetX() [1/2]

void PndFTSCATrackParamVector::SetX ( const float_v &  v)
inline

Definition at line 1416 of file PndFTSCATrackParamVector.h.

References v.

Referenced by PndFTSCAGBTracker::Refit_1(), and Rotate().

1416 { fX = v; }
__m128 v
Definition: P4_F32vec4.h:15

◆ SetX() [2/2]

void PndFTSCATrackParamVector::SetX ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1419 of file PndFTSCATrackParamVector.h.

References v.

1419 { fX(m) = v; }
__m128 m
Definition: P4_F32vec4.h:38
__m128 v
Definition: P4_F32vec4.h:15

◆ SetY() [1/2]

void PndFTSCATrackParamVector::SetY ( const float_v &  v)
inline

Definition at line 1417 of file PndFTSCATrackParamVector.h.

References v.

Referenced by PndFTSCAGBTracker::Refit_1(), and Rotate().

1417 { fP[0] = v; }
__m128 v
Definition: P4_F32vec4.h:15

◆ SetY() [2/2]

void PndFTSCATrackParamVector::SetY ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1420 of file PndFTSCATrackParamVector.h.

References m.

1420 { fP[0](m) = v; }
__m128 m
Definition: P4_F32vec4.h:38
__m128 v
Definition: P4_F32vec4.h:15

◆ SetZ() [1/2]

void PndFTSCATrackParamVector::SetZ ( const float_v &  v)
inline

Definition at line 1418 of file PndFTSCATrackParamVector.h.

References v.

1418 { fP[1] = v; }
__m128 v
Definition: P4_F32vec4.h:15

◆ SetZ() [2/2]

void PndFTSCATrackParamVector::SetZ ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1421 of file PndFTSCATrackParamVector.h.

References m.

1421 { fP[1](m) = v; }
__m128 m
Definition: P4_F32vec4.h:38
__m128 v
Definition: P4_F32vec4.h:15

◆ SignCosPhi()

float_v PndFTSCATrackParamVector::SignCosPhi ( ) const
inline

The sign of cos phi is always positive in the slice tracker. Only after coordinate transformation can the sign change to negative.

Definition at line 1344 of file PndFTSCATrackParamVector.h.

Referenced by PndFTSCATrackLinearisationVector::PndFTSCATrackLinearisationVector(), SetTrackParam(), SetTrackParamOne(), Tx1(), and Tx2().

1344 { return fSignCosPhi; }

◆ SinPhi()

float_v PndFTSCATrackParamVector::SinPhi ( ) const
inline

Definition at line 1327 of file PndFTSCATrackParamVector.h.

Referenced by GetCosPhi(), GetCosPhiPositive(), Rotate(), RotateXY(), TransportJ0ToX0(), TransportToX0(), Tx1(), and Tx2().

1327 { return fP[2]; }

◆ Transport() [1/3]

float_m PndFTSCATrackParamVector::Transport ( const int_v &  ista,
const PndFTSCAParam param,
const float_m &  mask = float_m(true) 
)

Referenced by SetErr2QPt().

◆ Transport() [2/3]

float_m PndFTSCATrackParamVector::Transport ( const FTSCAHitV hit,
const PndFTSCAParam p,
const float_m &  mask = float_m(true) 
)

◆ Transport() [3/3]

float_m PndFTSCATrackParamVector::Transport ( const FTSCAHit hit,
const PndFTSCAParam p,
const float_m &  mask = float_m(true) 
)

◆ TransportJ0ToX0()

float_m PndFTSCATrackParamVector::TransportJ0ToX0 ( const float_v &  x0,
const float_v &  cBz,
float_v &  extrDy,
float_v &  extrDz,
float_v  J[],
const float_m &  active = float_m(true) 
)
inline

Definition at line 2001 of file PndFTSCATrackParamVector.h.

References CAMath::Abs(), DzDs(), operator<<, operator>>, CAMath::RSqrt(), SinPhi(), and X().

Referenced by SetErr2QPt().

2002 {
2003  const float_v &ey = SinPhi();
2004  const float_v &dx = x - X();
2005  const float_v &exi = CAMath::RSqrt(float_v(Vc::One) - ey * ey); // RSqrt
2006 
2007  const float_v &dxBz = dx * cBz;
2008  const float_v &dS = dx * exi;
2009  const float_v &h2 = dS * exi * exi;
2010  const float_v &h4 = .5f * h2 * dxBz;
2011 
2012  float_m mask = active && CAMath::Abs(exi) <= 1.e4f;
2013 
2014  extrDy(active) = dS * ey;
2015  extrDz(active) = dS * DzDs();
2016  J[0](active) = dS;
2017  J[1](active) = 0;
2018  J[2](active) = h4;
2019  J[3](active) = dS;
2020  J[4](active) = dS;
2021  J[5](active) = 0;
2022  return active;
2023 }
static T Abs(const T &x)
Definition: PndCAMath.h:68
static T RSqrt(const T &x)
Definition: PndCAMath.h:62

◆ TransportToX0() [1/3]

float_m PndFTSCATrackParamVector::TransportToX0 ( const float_v &  x,
const float_v &  Bz,
const float  maxSinPhi = .999f,
const float_m &  mask = float_m(true) 
)

Referenced by SetErr2QPt().

◆ TransportToX0() [2/3]

float_m PndFTSCATrackParamVector::TransportToX0 ( const float_v &  x,
PndFTSCATrackLinearisationVector t0,
const float_v &  Bz,
const float  maxSinPhi = .999f,
float_v *  DL = 0,
const float_m &  mask = float_m(true) 
)

◆ TransportToX0() [3/3]

float_m PndFTSCATrackParamVector::TransportToX0 ( const float_v &  x,
const float_v &  sinPhi0,
const float_v &  Bz,
const float_v  maxSinPhi = .999f,
const float_m &  mask = float_m(true) 
)
inline

mvz start 23.01.2010

mvz end 23.01.2010

Definition at line 1518 of file PndFTSCATrackParamVector.h.

References CAMath::Abs(), debugKF(), DzDs(), QPt(), CAMath::RSqrt(), SinPhi(), and X().

1519 {
1520  //* Transport the track parameters to X=x, using linearization at phi0 with 0 curvature,
1521  //* and the field value Bz
1522  //* maxSinPhi is the max. allowed value for |t0.SinPhi()|
1523  //* linearisation of trajectory t0 is also transported to X=x,
1524  //* returns 1 if OK
1525  //*
1526 
1527  debugKF() << "Start TransportToX0(" << x << ", " << _mask << ")\n" << *this << std::endl;
1528 
1529  const float_v &ey = sinPhi0;
1530  const float_v &dx = x - X();
1531  const float_v &exi = float_v(Vc::One) * CAMath::RSqrt(float_v(Vc::One) - ey * ey); // RSqrt
1532 
1533  const float_v &dxBz = dx * Bz;
1534  const float_v &dS = dx * exi;
1535  const float_v &h2 = dS * exi * exi;
1536  const float_v &h4 = .5f * h2 * dxBz;
1537 //#define LOSE_DEBUG
1538 #ifdef LOSE_DEBUG
1539  std::cout << " TrTo-sinPhi0 = " << sinPhi0 << std::endl;
1540 #endif
1541  // const float_v &sinPhi = SinPhi() * (float_v( Vc::One ) - 0.5f * dxBz * QPt() *dxBz * QPt()/ ( float_v( Vc::One ) - SinPhi()*SinPhi() )) + dxBz * QPt();
1543  const float_v &sinPhi = SinPhi() + dxBz * QPt();
1545 #ifdef LOSE_DEBUG
1546  std::cout << " TrTo-sinPhi = " << sinPhi << std::endl;
1547 #endif
1548  float_m mask = _mask && CAMath::Abs(exi) <= 1.e4f;
1549  mask &= ((CAMath::Abs(sinPhi) <= maxSinPhi) || (maxSinPhi <= 0.f));
1550 
1551  fX(mask) += dx;
1552  fP[0](mask) += dS * ey + h2 * (SinPhi() - ey) + h4 * QPt();
1553  fP[1](mask) += dS * DzDs();
1554  fP[2](mask) = sinPhi;
1555 
1556  // const float_v c00 = fC[0];
1557  // const float_v c11 = fC[2];
1558  const float_v c20 = fC[3];
1559  // const float_v c21 = fC[4];
1560  const float_v c22 = fC[5];
1561  // const float_v c30 = fC[6];
1562  const float_v c31 = fC[7];
1563  // const float_v c32 = fC[8];
1564  const float_v c33 = fC[9];
1565  const float_v c40 = fC[10];
1566  // const float_v c41 = fC[11];
1567  const float_v c42 = fC[12];
1568  // const float_v c43 = fC[13];
1569  const float_v c44 = fC[14];
1570 
1571  const float_v two(2.f);
1572 
1573  fC[0](mask) += h2 * h2 * c22 + h4 * h4 * c44 + two * (h2 * c20 + h4 * (c40 + h2 * c42));
1574 
1575  // fC[1] ( mask ) += h2 * c21 + h4 * c41 + dS * ( c30 + h2 * c32 + h4 * c43 );
1576  fC[2](mask) += dS * (two * c31 + dS * c33);
1577 
1578  fC[3](mask) += h2 * c22 + h4 * c42 + dxBz * (c40 + h2 * c42 + h4 * c44);
1579  // fC[4] ( mask ) += dS * c32 + dxBz * ( c41 + dS * c43 );
1580  const float_v &dxBz_c44 = dxBz * c44;
1581  fC[12](mask) += dxBz_c44;
1582  fC[5](mask) += dxBz * (two * c42 + dxBz_c44);
1583 
1584  // fC[6] ( mask ) += h2 * c32 + h4 * c43;
1585  fC[7](mask) += dS * c33;
1586  // fC[8] ( mask ) += dxBz * c43;
1587  // fC[9] ( mask ) = c33;
1588 
1589  fC[10](mask) += h2 * c42 + h4 * c44;
1590  // fC[11]( mask ) += dS * c43;
1591  // fC[13]( mask ) = c43;
1592  // fC[14]( mask ) = c44;
1593 
1594  debugKF() << mask << "\n" << *this << std::endl;
1595  return mask;
1596 }
static T Abs(const T &x)
Definition: PndCAMath.h:68
static T RSqrt(const T &x)
Definition: PndCAMath.h:62
float f
Definition: P4_F32vec4.h:32
PndFTSCANoDebugStream & debugKF()
Definition: debug.h:138

◆ TransportToX0WithMaterial() [1/3]

float_m PndFTSCATrackParamVector::TransportToX0WithMaterial ( const float_v &  x,
const float_v &  XOverX0,
const float_v &  XThimesRho,
const float_v &  Bz,
const float  maxSinPhi = .999f 
)

Referenced by SetErr2QPt().

◆ TransportToX0WithMaterial() [2/3]

float_m PndFTSCATrackParamVector::TransportToX0WithMaterial ( const float_v &  x,
PndFTSCATrackLinearisationVector t0,
PndFTSCATrackFitParam par,
const float_v &  XOverX0,
const float_v &  XThimesRho,
const float_v &  Bz,
const float  maxSinPhi = .999f,
const float_m &  mask = float_m(true) 
)

◆ TransportToX0WithMaterial() [3/3]

float_m PndFTSCATrackParamVector::TransportToX0WithMaterial ( const float_v &  x,
PndFTSCATrackFitParam par,
const float_v &  XOverX0,
const float_v &  XThimesRho,
const float_v &  Bz,
const float  maxSinPhi = .999f 
)

◆ Tx1()

float_v PndFTSCATrackParamVector::Tx1 ( ) const
inline

Definition at line 1336 of file PndFTSCATrackParamVector.h.

References SignCosPhi(), SinPhi(), and sqrt().

1336 { return SinPhi() / (SignCosPhi() * sqrt(1 - SinPhi() * SinPhi())); } // CHECKME
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:40

◆ Tx2()

float_v PndFTSCATrackParamVector::Tx2 ( ) const
inline

Definition at line 1337 of file PndFTSCATrackParamVector.h.

References DzDs(), SignCosPhi(), SinPhi(), and sqrt().

1337 { return DzDs() / (SignCosPhi() * sqrt(1 - SinPhi() * SinPhi())); } // dx2/dx0 = dz/dx
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:40

◆ X()

float_v PndFTSCATrackParamVector::X ( ) const
inline

◆ X0()

float_v PndFTSCATrackParamVector::X0 ( ) const
inline

Definition at line 1333 of file PndFTSCATrackParamVector.h.

References X().

1333 { return X(); }

◆ X1()

float_v PndFTSCATrackParamVector::X1 ( ) const
inline

Definition at line 1334 of file PndFTSCATrackParamVector.h.

References Y().

1334 { return Y(); }

◆ X2()

float_v PndFTSCATrackParamVector::X2 ( ) const
inline

Definition at line 1335 of file PndFTSCATrackParamVector.h.

References Z().

1335 { return Z(); }

◆ Y()

float_v PndFTSCATrackParamVector::Y ( ) const
inline

Definition at line 1325 of file PndFTSCATrackParamVector.h.

Referenced by FilterVtx(), PndFTSCAGBTracker::Refit_1(), Rotate(), RotateXY(), and X1().

1325 { return fP[0]; }

◆ Z()

float_v PndFTSCATrackParamVector::Z ( ) const
inline

Definition at line 1326 of file PndFTSCATrackParamVector.h.

Referenced by FilterVtx(), PndFTSCAGBTracker::Refit_1(), and X2().

1326 { return fP[1]; }

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  ,
const PndFTSCATrackParamVector  
)
friend

Referenced by TransportJ0ToX0().

◆ operator>>

std::istream& operator>> ( std::istream &  ,
PndFTSCATrackParamVector  
)
friend

Referenced by TransportJ0ToX0().

◆ PndFTSCATrackParam

friend class PndFTSCATrackParam
friend

Definition at line 1378 of file PndFTSCATrackParamVector.h.


The documentation for this class was generated from the following file: