FairRoot/PandaRoot
Classes | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
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)
 

Private Member Functions

const float_v * Par () const
 
const float_v * Cov () const
 
float_v * Par ()
 
float_v * Cov ()
 

Private Attributes

float_v fX
 
float_v fSignCosPhi
 
float_v fP [5]
 
float_v fC [15]
 
float_v fChi2
 
int_v fNDF
 
float_v fAlpha
 

Friends

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

Detailed Description

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

Definition at line 1269 of file PndFTSCATrackParamVector.h.

Constructor & Destructor Documentation

PndFTSCATrackParamVector::PndFTSCATrackParamVector ( )
inline

Definition at line 1274 of file PndFTSCATrackParamVector.h.

References fP, and i.

1275  : fX( Vc::Zero ),
1276  fSignCosPhi( Vc::Zero ),
1277  fChi2( Vc::Zero ),
1278  fNDF( Vc::Zero )
1279  {
1280  for ( int i = 0; i < 5; ++i ) fP[i].setZero();
1281  for ( int i = 0; i < 15; ++i ) fC[i].setZero();
1282  }
Int_t i
Definition: run_full.C:25
static const fvec Zero

Member Function Documentation

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

Definition at line 1504 of file PndFTSCATrackParamVector.cxx.

References FTSCATarget::B(), FTSCATarget::Err2X1(), FTSCATarget::Err2X2(), FilterVtx(), fNDF, FTSCATarget::NDF(), TransportJ0ToX0(), FTSCATarget::X0(), X0(), FTSCATarget::X1(), and FTSCATarget::X2().

Referenced by PndFTSCAGBTracker::FitTrackCA().

1505 {
1506  float_m active = mask;
1507 
1508  float_v eY, eZ;
1509  float_v dz = target.X0() - X0();
1510  float_v J[6];
1511  // H = 1 0 J[0] J[1] J[2]
1512  // 0 1 J[3] J[4] J[5]
1513  active &= TransportJ0ToX0( target.X0(), static_cast<float_v>(target.B()), eY, eZ, J, active );
1514  active &= FilterVtx( target.X1(), target.X2(), CAX1X2MeasurementInfo( target.Err2X1(), 0.f, target.Err2X2() ), eY, eZ, J, active );
1515  fNDF (static_cast<int_m>(active)) += target.NDF();
1516 
1517  return active;
1518 }
float X1() const
Definition: FTSCATarget.h:39
int NDF() const
Definition: FTSCATarget.h:48
float Err2X1() const
Definition: FTSCATarget.h:42
const CAFieldValue & B(int i=0) const
Definition: FTSCATarget.h:46
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 X2() const
Definition: FTSCATarget.h:40
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 X0() const
Definition: FTSCATarget.h:38
float Err2X2() const
Definition: FTSCATarget.h:43
float_v PndFTSCATrackParamVector::Angle ( ) const
inline
float_v PndFTSCATrackParamVector::ApproximateBetheBloch ( const float_v &  beta2)
static

Definition at line 1019 of file PndFTSCATrackParamVector.cxx.

References CAMath::Log(), and One.

Referenced by CalculateFitParameters().

1020 {
1021  //------------------------------------------------------------------
1022  // This is an approximation of the Bethe-Bloch formula with
1023  // the density effect taken into account at beta*gamma > 3.5
1024  // (the approximation is reasonable only for solid materials)
1025  //------------------------------------------------------------------
1026 
1027  const float_v &beta2_1subBeta2 = beta2 / ( float_v( Vc::One ) - beta2 ); // beta2 * CAMath::Reciprocal( float_v( Vc::One ) - beta2 );
1028  const float_v &_0p000153_beta2 = 0.153e-3f / beta2;
1029  const float log_3p5mul5940 = 9.942227380852058f; // log( 3.5 * 5940 )
1030  const float log_5940 = 8.68946441235669f; // log( 5940 )
1031  const float_v log_beta2_1subBeta2 = CAMath::Log( beta2_1subBeta2 );
1032 
1033  float_v ret = _0p000153_beta2 * ( log_5940 + log_beta2_1subBeta2 - beta2 );
1034  ret( beta2_1subBeta2 > 3.5f*3.5f ) =
1035  _0p000153_beta2 * ( log_3p5mul5940 + 0.5f * log_beta2_1subBeta2 - beta2 );
1036  ret.setZero( beta2 >= float_v( Vc::One ) );
1037  return ret;
1038 }
TFile * f
Definition: bump_analys.C:12
static T Log(const T &x)
Definition: PndCAMath.h:40
static const fvec One
float_v PndFTSCATrackParamVector::BetheBlochGas ( const float_v &  bg)
static

Definition at line 998 of file PndFTSCATrackParamVector.cxx.

References BetheBlochGeant(), f, and x0.

999 {
1000  //------------------------------------------------------------------
1001  // This is an approximation of the Bethe-Bloch formula,
1002  // reasonable for gas materials.
1003  // All the parameters are, in fact, for Ne.
1004  // The returned value is in [GeV]
1005  //------------------------------------------------------------------
1006 
1007  const float_v rho = 0.9e-3f;
1008  const float_v x0 = 2.f;
1009  const float_v x1 = 4.f;
1010  const float_v mI = 140.e-9f;
1011  const float_v mZA = 0.49555f;
1012 
1013  return BetheBlochGeant( bg, rho, x0, x1, mI, mZA );
1014 }
Double_t x0
Definition: checkhelixhit.C:70
TFile * f
Definition: bump_analys.C:12
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)
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

Definition at line 945 of file PndFTSCATrackParamVector.cxx.

References f, CAMath::Log(), mK, One, r, CAMath::Sqrt(), x, x0, and Zero.

Referenced by BetheBlochGas(), and BetheBlochSolid().

951 {
952  //
953  // This is the parameterization of the Bethe-Bloch formula inspired by Geant.
954  //
955  // bg2 - (beta*gamma)^2
956  // kp0 - density [g/cm^3]
957  // kp1 - density effect first junction point
958  // kp2 - density effect second junction point
959  // kp3 - mean excitation energy [GeV]
960  // kp4 - mean Z/A
961  //
962  // The default values for the kp* parameters are for silicon.
963  // The returned value is in [GeV/(g/cm^2)].
964  //
965 
966  const float mK = 0.307075e-3f; // [GeV*cm^2/g]
967  const float _2me = 1.022e-3f; // [GeV/c^2]
968  const float_v &rho = kp0;
969  const float_v &x0 = kp1 * 2.303f;
970  const float_v &x1 = kp2 * 2.303f;
971  const float_v &mI = kp3;
972  const float_v &mZA = kp4;
973  const float_v &maxT = _2me * bg2; // neglecting the electron mass
974 
975  //*** Density effect
976  float_v d2( Vc::Zero );
977  const float_v x = 0.5f * CAMath::Log( bg2 );
978  const float_v lhwI = CAMath::Log( 28.816f * 1e-9f * CAMath::Sqrt( rho * mZA ) / mI );
979  d2( x > x1 ) = lhwI + x - 0.5f;
980  const float_v &r = ( x1 - x ) / ( x1 - x0 );
981  d2( x > x0 && x <= x1 ) = lhwI + x - 0.5f + ( 0.5f - lhwI - x0 ) * r * r * r;
982 
983  return mK*mZA*( float_v( Vc::One ) + bg2 ) / bg2*( 0.5f*CAMath::Log( maxT * maxT / ( mI*mI ) ) - bg2 / ( float_v( Vc::One ) + bg2 ) - d2 );
984 }
Double_t x0
Definition: checkhelixhit.C:70
double mK
double r
Definition: RiemannTest.C:14
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
static const fvec Zero
TFile * f
Definition: bump_analys.C:12
Double_t x
static T Log(const T &x)
Definition: PndCAMath.h:40
static const fvec One
float_v PndFTSCATrackParamVector::BetheBlochSolid ( const float_v &  bg)
static

Definition at line 986 of file PndFTSCATrackParamVector.cxx.

References BetheBlochGeant().

987 {
988  //------------------------------------------------------------------
989  // This is an approximation of the Bethe-Bloch formula,
990  // reasonable for solid materials.
991  // All the parameters are, in fact, for Si.
992  // The returned value is in [GeV]
993  //------------------------------------------------------------------
994 
995  return BetheBlochGeant( bg );
996 }
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)
void PndFTSCATrackParamVector::CalculateFitParameters ( PndFTSCATrackFitParam par,
const float_v &  mass = 0.13957f 
)

Definition at line 1041 of file PndFTSCATrackParamVector.cxx.

References ApproximateBetheBloch(), PndFTSCATrackParamVector::PndFTSCATrackFitParam::fBethe, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fE, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fEP2, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fK22, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fK33, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fK43, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fK44, fP, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fSigmadE2, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fTheta2, One, p2, and CAMath::Sqrt().

Referenced by PndFTSCAGBTracker::FitTrack(), Transport(), and TransportToX0WithMaterial().

1042 {
1043  //*!
1044 
1045  const float_v p2 = ( float_v( Vc::One ) + fP[3] * fP[3] );
1046  const float_v k2 = fP[4] * fP[4];
1047  const float_v mass2 = mass * mass;
1048  const float_v beta2 = p2 / ( p2 + mass2 * k2 );
1049 
1050  float_v pp2 = 10000.f; pp2( k2 > 1.e-8f ) = p2 / k2; // impuls 2
1051 
1052  //par.fBethe = BetheBlochGas( pp2/mass2);
1053  // par.fBethe = ApproximateBetheBloch( pp2 / mass2 ); Why pp2 / mass2?, should be beta2.
1054  par.fBethe = ApproximateBetheBloch( beta2 );
1055  par.fE = CAMath::Sqrt( pp2 + mass2 );
1056  par.fTheta2 = 14.1f * 14.1f / ( beta2 * pp2 * 1.e6f );
1057  par.fEP2 = par.fE / pp2;
1058 
1059  // Approximate energy loss fluctuation (M.Ivanov)
1060 
1061  const float knst = 0.07f; // To be tuned.
1062  par.fSigmadE2 = knst * par.fEP2 * fP[4];
1063  par.fSigmadE2 = par.fSigmadE2 * par.fSigmadE2;
1064 
1065  par.fK22 = ( float_v( Vc::One ) + fP[3] * fP[3] );
1066  par.fK33 = par.fK22 * par.fK22;
1067  par.fK43 = fP[3] * fP[4] * par.fK22;
1068  par.fK44 = fP[3] * fP[3] * fP[4] * fP[4];
1069 }
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
Double_t par[3]
static float_v ApproximateBetheBloch(const float_v &beta2)
TFile * f
Definition: bump_analys.C:12
TPad * p2
Definition: hist-t7.C:117
static const fvec One
float_v PndFTSCATrackParamVector::Chi2 ( ) const
inline

Definition at line 1331 of file PndFTSCATrackParamVector.h.

References fChi2.

Referenced by ConvertTrackParamToVector(), operator<<(), and PndFTSCAGBTracker::Refit_1().

1331 { return fChi2; }
void PndFTSCATrackParamVector::ConvertTrackParamToVector ( PndFTSCATrackParam  t0[float_v::Size],
int  nTracksV 
)

Definition at line 557 of file PndFTSCATrackParamVector.cxx.

References Angle(), Chi2(), Cov(), NDF(), Par(), SetAngle(), SetChi2(), SetCov(), SetNDF(), SetPar(), SetSignCosPhi(), SetX(), SignCosPhi(), t0, and X().

Referenced by PndFTSCAGBTracker::FitTracks().

559 {
560  float_v tmpVec;
561  int_v tmpVecShort;
562  float_v::Memory tmpFloat;
563  int_v::Memory tmpShort;
564 
565  for(int iV=0; iV < nTracksV; iV++) tmpFloat[iV] = t0[iV].X();
566  tmpVec.load( tmpFloat );
567  SetX(tmpVec);
568  for(int iV=0; iV < nTracksV; iV++) tmpFloat[iV] = t0[iV].SignCosPhi();
569  tmpVec.load( tmpFloat );
570  SetSignCosPhi(tmpVec);
571 
572  for(int iP=0; iP<5; iP++)
573  {
574  for(int iV=0; iV < nTracksV; iV++) tmpFloat[iV] = t0[iV].Par()[iP];
575  tmpVec.load( tmpFloat );
576  SetPar(iP,tmpVec);
577  }
578  for(int iC=0; iC<15; iC++)
579  {
580  for(int iV=0; iV < nTracksV; iV++) tmpFloat[iV] = t0[iV].Cov()[iC];
581  tmpVec.load( tmpFloat );
582  SetCov(iC,tmpVec);
583  }
584  for(int iV=0; iV < nTracksV; iV++) tmpFloat[iV] = t0[iV].Chi2();
585  tmpVec.load( tmpFloat );
586  SetChi2(tmpVec);
587  for(int iV=0; iV < nTracksV; iV++) tmpShort[iV] = t0[iV].NDF();
588  tmpVecShort.load( tmpShort );
589  SetNDF(tmpVecShort);
590 
591  for(int iV=0; iV < nTracksV; iV++) tmpFloat[iV] = t0[iV].Angle();
592  tmpVec.load( tmpFloat );
593  SetAngle(tmpVec);
594 }
void SetCov(int i, const float_v &v)
const float_v * Par() const
void SetSignCosPhi(const float_v &v)
const float_v * Cov() const
void SetPar(int i, const float_v &v)
float_m PndFTSCATrackParamVector::CorrectForMeanMaterial ( const float_v &  xOverX0,
const float_v &  xTimesRho,
const PndFTSCATrackFitParam par,
const float_m &  _mask 
)

Definition at line 1072 of file PndFTSCATrackParamVector.cxx.

References CAMath::Abs(), corr, dE, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fBethe, fC, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fE, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fEP2, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fK22, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fK33, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fK43, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fK44, fP, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fSigmadE2, PndFTSCATrackParamVector::PndFTSCATrackFitParam::fTheta2, One, and Zero.

Referenced by TransportToX0WithMaterial().

1073 {
1074  //------------------------------------------------------------------
1075  // This function corrects the track parameters for the crossed material.
1076  // "xOverX0" - X/X0, the thickness in units of the radiation length.
1077  // "xTimesRho" - is the product length*density (g/cm^2).
1078  //------------------------------------------------------------------
1079 
1080  float_v &fC22 = fC[5];
1081  float_v &fC33 = fC[9];
1082  float_v &fC40 = fC[10];
1083  float_v &fC41 = fC[11];
1084  float_v &fC42 = fC[12];
1085  float_v &fC43 = fC[13];
1086  float_v &fC44 = fC[14];
1087 
1088  //Energy losses************************
1089 
1090  const float_v &dE = par.fBethe * xTimesRho;
1091 // assert( (dE >= float_v(Vc::Zero) || !_mask ).isFull() );
1092  float_m mask = _mask && dE <= 0.3f * par.fE; //30% energy loss is too much!
1093  const float_v &corr = ( float_v( Vc::One ) - par.fEP2 * dE );
1094  mask &= corr >= 0.3f && corr <= 1.3f;
1095 
1096  fP[4]( mask ) *= corr;
1097  fC40 ( mask ) *= corr;
1098  fC41 ( mask ) *= corr;
1099  fC42 ( mask ) *= corr;
1100  fC43 ( mask ) *= corr;
1101  fC44 ( mask ) *= corr * corr;
1102  fC44 ( mask ) += par.fSigmadE2 * CAMath::Abs( dE );
1103 
1104  //Multiple scattering******************
1105 
1106  assert( (xOverX0 >= float_v(Vc::Zero) || !mask ).isFull() );
1107  const float_v &theta2 = par.fTheta2 * xOverX0;
1108  fC22( mask ) += theta2 * par.fK22 * ( float_v( Vc::One ) - fP[2] * fP[2] );
1109  fC33( mask ) += theta2 * par.fK33;
1110  fC43( mask ) += theta2 * par.fK43;
1111  fC44( mask ) += theta2 * par.fK44;
1112 
1113  return mask;
1114 }
PndPidCorrelator * corr
Double_t par[3]
static const fvec Zero
static T Abs(const T &x)
Definition: PndCAMath.h:39
Double_t dE
Definition: anasim.C:58
static const fvec One
const float_v& PndFTSCATrackParamVector::Cov ( int  i) const
inline
const float_v* PndFTSCATrackParamVector::Cov ( ) const
inlineprivate

Definition at line 1366 of file PndFTSCATrackParamVector.h.

References fC.

Referenced by ConvertTrackParamToVector().

1366 { return fC; }
float_v* PndFTSCATrackParamVector::Cov ( )
inlineprivate

Definition at line 1368 of file PndFTSCATrackParamVector.h.

References fC.

1368 { return fC; }
float_v PndFTSCATrackParamVector::DzDs ( ) const
inline

Definition at line 1314 of file PndFTSCATrackParamVector.h.

References fP.

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

1314 { return fP[3]; }
float_v PndFTSCATrackParamVector::Err2DzDs ( ) const
inline

Definition at line 1337 of file PndFTSCATrackParamVector.h.

1337 { return fC[9]; }
float_v PndFTSCATrackParamVector::Err2QMomentum ( ) const
inline

Definition at line 1358 of file PndFTSCATrackParamVector.h.

1358 { return fC[14]; }
float_v PndFTSCATrackParamVector::Err2QPt ( ) const
inline

Definition at line 1338 of file PndFTSCATrackParamVector.h.

1338 { return fC[14]; }
float_v PndFTSCATrackParamVector::Err2SinPhi ( ) const
inline

Definition at line 1336 of file PndFTSCATrackParamVector.h.

1336 { return fC[5]; }
float_v PndFTSCATrackParamVector::Err2X1 ( ) const
inline

Definition at line 1354 of file PndFTSCATrackParamVector.h.

Referenced by PndFTSCAGBTracker::PickUpHits().

1354 { return fC[0]; }
float_v PndFTSCATrackParamVector::Err2X2 ( ) const
inline

Definition at line 1355 of file PndFTSCATrackParamVector.h.

1355 { return fC[2]; }
float_v PndFTSCATrackParamVector::Err2Y ( ) const
inline

Definition at line 1334 of file PndFTSCATrackParamVector.h.

1334 { return fC[0]; }
float_v PndFTSCATrackParamVector::Err2Z ( ) const
inline

Definition at line 1335 of file PndFTSCATrackParamVector.h.

1335 { return fC[2]; }
float_m PndFTSCATrackParamVector::Filter ( const FTSCAHitV hit,
const PndFTSCAParam param,
const float_m &  mask = float_m( true ),
const float_v &  chi2Cut = 10e10f 
)

Definition at line 1464 of file PndFTSCATrackParamVector.cxx.

References FTSCAHitV::Err2X1(), FTSCAHitV::Err2X2(), FTSCAHitV::ErrX12(), FTSCAStation::f, FilterWithMaterial(), FTSCAHitV::IStations(), FTSCAHitV::IsValid(), FTSCAStation::NDF, NDF(), PndFTSCAParam::Station(), FTSCAHitV::X1(), and FTSCAHitV::X2().

Referenced by CAFunctionality::FitIteration(), PndFTSCAGBTracker::FitTrackCA(), PndFTSCAGBTracker::PickUpHits(), and PndFTSCAGBTracker::Refit().

1465 {
1466  const int iSta = hit.IStations()[hit.IsValid().firstOne()];
1467  if ( param.Station( iSta ).NDF == 1 ) {
1468 #ifdef DRIFT_TUBES
1469  return FilterWithMaterial( hit.XWire1(), hit.XWire2(), hit.RSigned(), param.Station( iSta ).f, hit.Err2R(), 0.999f, mask, chi2Cut );
1470 #else
1471  assert(0);
1472 #endif
1473  }
1474  return FilterWithMaterial( hit.X1(), hit.X2(), hit.Err2X1(), hit.ErrX12(), hit.Err2X2(), 0.999f, mask, int_v(param.Station( hit.IStations()[0] ).NDF), chi2Cut );
1475 }
float_v X2() const
Definition: FTSCAHitsV.h:51
float_v Err2X1() const
Definition: FTSCAHitsV.h:56
const FTSCAStation & Station(short i) const
Definition: PndFTSCAParam.h:45
FTSCAStripInfo f
Definition: FTSCAStation.h:32
float_m IsValid() const
Definition: FTSCAHitsV.h:38
float_v X1() const
Definition: FTSCAHitsV.h:50
int_v IStations() const
Definition: FTSCAHitsV.h:45
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_v Err2X2() const
Definition: FTSCAHitsV.h:58
float_v ErrX12() const
Definition: FTSCAHitsV.h:57
float_m PndFTSCATrackParamVector::Filter ( const FTSCAHit hit,
const PndFTSCAParam param,
const float_m &  mask = float_m( true ),
const float_v &  chi2Cut = 10e10f 
)

Definition at line 1492 of file PndFTSCATrackParamVector.cxx.

References FTSCAHit::Err2X1(), FTSCAHit::Err2X2(), FTSCAHit::ErrX12(), FilterWithMaterial(), FTSCAHit::IStation(), NDF(), PndFTSCAParam::Station(), FTSCAHit::X1(), and FTSCAHit::X2().

1493 {
1494  if ( param.Station( hit.IStation() ).NDF == 1 ) {
1495 #ifdef DRIFT_TUBES
1496  return FilterWithMaterial( hit.XWire1(), hit.XWire2(), hit.RSigned(), param.Station( hit.IStation() ).f, hit.Err2R(), 0.999f, mask, chi2Cut );
1497 #else
1498  assert(0);
1499 #endif
1500  }
1501  return FilterWithMaterial( hit.X1(), hit.X2(), hit.Err2X1(), hit.ErrX12(), hit.Err2X2(), 0.999f, mask, int_v(param.Station( hit.IStation() ).NDF), chi2Cut );
1502 }
float X2() const
Definition: FTSCAHits.h:43
const FTSCAStation & Station(short i) const
Definition: PndFTSCAParam.h:45
char IStation() const
Definition: FTSCAHits.h:33
float X1() const
Definition: FTSCAHits.h:42
float Err2X2() const
Definition: FTSCAHits.h:52
TFile * f
Definition: bump_analys.C:12
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 Err2X1() const
Definition: FTSCAHits.h:50
float ErrX12() const
Definition: FTSCAHits.h:51
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 1917 of file PndFTSCATrackParamVector.h.

References CAX1X2MeasurementInfo::C00(), CAX1X2MeasurementInfo::C10(), c10, CAX1X2MeasurementInfo::C11(), c11, c20, c21, c22, fChi2, fP, Y(), and Z().

Referenced by AddTarget().

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

Definition at line 1168 of file PndFTSCATrackParamVector.cxx.

References CAMath::Abs(), c10, c11, c20, c21, d, f, fC, fChi2, fNDF, fP, One, and mrfdata_8b_error::success.

Referenced by Filter(), and PndFTSCAGBTracker::FitTrack().

1169 {
1170  assert( maxSinPhi > 0.f );
1171  //* Add the y,z measurement with the Kalman filter
1172 
1173  const float_v c00 = fC[0];
1174  const float_v c10 = fC[1];
1175  const float_v c11 = fC[2];
1176  const float_v c20 = fC[3];
1177  const float_v c21 = fC[4];
1178 // float c22 = fC[5];
1179  const float_v c30 = fC[6];
1180  const float_v c31 = fC[7];
1181 // float c32 = fC[8];
1182 // float c33 = fC[9];
1183  const float_v c40 = fC[10];
1184  const float_v c41 = fC[11];
1185 // float c42 = fC[12];
1186 // float c43 = fC[13];
1187 // float c44 = fC[14];
1188 
1189 
1190  const float_v
1191  z0 = y - fP[0],
1192  z1 = z - fP[1];
1193 
1194  // foreach_bit( int i, mask ) {
1195  // cout << i << " e = " << err2Y[i] << " " << errYZ[i] << " " << err2Z[i] << " c = " << c00[i] << " " << c10[i] << " " << c11[i] << " z = " << z0[i] << " " << z1[i] << " chi2 = " << fChi2[i]<< endl;
1196  // }
1197 
1198  err2Y += c00;
1199  err2Z += c11;
1200  errYZ += c10;
1201  float_v d = float_v( Vc::One ) / ( err2Y*err2Z - errYZ*errYZ );
1202 
1203 
1204  float_m success = mask;// if ( ISUNLIKELY( err2Y < 1.e-8f ) || ISUNLIKELY( err2Z < 1.e-8f ) ) return 0;
1205  success &= (err2Y > 1.e-8f ) && ( err2Z > 1.e-8f );
1206 
1207  const float_v mS0 = err2Z*d;
1208  const float_v mS1 = -errYZ*d;
1209  const float_v mS2 = err2Y*d;
1210 
1211  // K = CHtS
1212 
1213  const float_v
1214  k00 = c00 * mS0 + c10*mS1, k01 = c00 * mS1 + c10*mS2,
1215  k10 = c10 * mS0 + c11*mS1, k11 = c10 * mS1 + c11*mS2,
1216  k20 = c20 * mS0 + c21*mS1, k21 = c20 * mS1 + c21*mS2,
1217  k30 = c30 * mS0 + c31*mS1, k31 = c30 * mS1 + c31*mS2,
1218  k40 = c40 * mS0 + c41*mS1, k41 = c40 * mS1 + c41*mS2;
1219 
1220  const float_v sinPhi = fP[2] + k20 * z0 + k21 * z1;
1221 
1222  success &= CAMath::Abs( sinPhi ) < maxSinPhi;
1223 
1224  fNDF( static_cast<int_m>(success) ) += hitNDF;
1225  fChi2(success) += mS0 * z0 * z0 + 2 * z0 * z1 * mS1 + mS2 * z1 * z1;
1226  success &= fChi2 < chi2Cut;
1227  if ( success.isEmpty() ) return success; // TODO move upper
1228 
1229  // foreach_bit( int i, success ) {
1230  // cout << i << " " << mS0[i] << " " << mS1[i] << " " << mS2[i] << " " << (mS0 * z0 * z0)[i] << " " << (2 * z0 * z1 * mS1)[i] << " " << (mS2 * z1 * z1)[i] << " " << fChi2[i] << endl;
1231  // }
1232  // std::cout << success << std::endl;
1233  fP[ 0](success) += k00 * z0 + k01 * z1;
1234  fP[ 1](success) += k10 * z0 + k11 * z1;
1235  fP[ 2](success) = sinPhi ;
1236  fP[ 3](success) += k30 * z0 + k31 * z1;
1237  fP[ 4](success) += k40 * z0 + k41 * z1;
1238  fC[ 0](success) -= (k00 * c00 + k01 * c10); //c00
1239 
1240  fC[ 1](success) -= (k10 * c00 + k11 * c10); //c10
1241  fC[ 2](success) -= (k10 * c10 + k11 * c11); //c11
1242 
1243  fC[ 3](success) -= (k20 * c00 + k21 * c10); //c20
1244  fC[ 4](success) -= (k20 * c10 + k21 * c11); //c21
1245  fC[ 5](success) -= (k20 * c20 + k21 * c21); //c22
1246 
1247  fC[ 6](success) -= (k30 * c00 + k31 * c10); //c30
1248  fC[ 7](success) -= (k30 * c10 + k31 * c11); //c31
1249  fC[ 8](success) -= (k30 * c20 + k31 * c21); //c32
1250  fC[ 9](success) -= (k30 * c30 + k31 * c31); //c33
1251 
1252  fC[10](success) -= (k40 * c00 + k41 * c10); //c40
1253  fC[11](success) -= (k40 * c10 + k41 * c11); //c41
1254 // std::cout << fC[11] << " " << fC[10] << std::endl;
1255 
1256  fC[12](success) -= (k40 * c20 + k41 * c21); //c42
1257  fC[13](success) -= (k40 * c30 + k41 * c31); //c43
1258  fC[14](success) -= (k40 * c40 + k41 * c41); //c44
1259 
1260  return success;
1261 }
Double_t z0
Definition: checkhelixhit.C:62
TCanvas * c11
TObjArray * d
TCanvas * c10
TCanvas * c21
static T Abs(const T &x)
Definition: PndCAMath.h:39
TFile * f
Definition: bump_analys.C:12
Double_t z
TCanvas * c20
Double_t y
static const fvec One
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 
)

Definition at line 1263 of file PndFTSCATrackParamVector.cxx.

References CAMath::Abs(), c10, c11, c20, c21, FTSCAStripInfo::cos, f, fC, fChi2, fNDF, fP, FTSCAStripInfo::sin, mrfdata_8b_error::success, and z.

1266 {
1267  assert( maxSinPhi > 0.f );
1268  float_m success = mask;
1269  success &= ( err2 > 1.e-8f );
1270 
1271  const float_v& c00 = fC[0];
1272  const float_v& c10 = fC[1];
1273  const float_v& c11 = fC[2];
1274 
1275  const float_v& sb = info.sin;
1276  const float_v& cb = info.cos;
1277 
1278  const float_v u = cb*y + sb*z;
1279  const float_v zeta = cb*fP[0] + sb*fP[1] - u;
1280 
1281  // F = CH'
1282  const float_v F0 = cb*c00 + sb*c10;
1283  const float_v F1 = cb*c10 + sb*c11;
1284 
1285  const float_v HCH = ( F0*cb + F1*sb );
1286 
1287  const float_v wi = 1/( err2 + HCH );
1288  const float_v zetawi = zeta *wi;
1289 
1290  fChi2(success) += zeta * zetawi;
1291  success &= fChi2 < chi2Cut;
1292  if ( success.isEmpty() ) return success;
1293 
1294  const float_v& c20 = fC[3];
1295  const float_v& c21 = fC[4];
1296  const float_v& c30 = fC[6];
1297  const float_v& c31 = fC[7];
1298  const float_v& c40 = fC[10];
1299  const float_v& c41 = fC[11];
1300 
1301  const float_v F2 = cb*c20 + sb*c21;
1302  const float_v F3 = cb*c30 + sb*c31;
1303  const float_v F4 = cb*c40 + sb*c41;
1304 
1305 
1306  const float_v K1 = F1*wi;
1307  const float_v K2 = F2*wi;
1308  const float_v K3 = F3*wi;
1309  const float_v K4 = F4*wi;
1310 
1311 
1312  const float_v sinPhi = fP[2] - F2*zetawi;
1313 
1314  success &= CAMath::Abs( sinPhi ) < maxSinPhi;
1315 
1316  fNDF( static_cast<int_m>(success) ) += 1; // TODO
1317 
1318  fP[ 0](success) -= F0*zetawi;
1319  fP[ 1](success) -= F1*zetawi;
1320  fP[ 2](success) = sinPhi ;
1321  fP[ 3](success) -= F3*zetawi;
1322  fP[ 4](success) -= F4*zetawi;
1323  fC[ 0](success) -= F0*F0*wi;
1324 
1325  fC[ 1](success) -= K1*F0;
1326  fC[ 2](success) -= K1*F1;
1327 
1328  fC[ 3](success) -= K2*F0;
1329  fC[ 4](success) -= K2*F1;
1330  fC[ 5](success) -= K2*F2;
1331 
1332  fC[ 6](success) -= K3*F0;
1333  fC[ 7](success) -= K3*F1;
1334  fC[ 8](success) -= K3*F2;
1335  fC[ 9](success) -= K3*F3;
1336 
1337  fC[10](success) -= K4*F0;
1338  fC[11](success) -= K4*F1;
1339 
1340  fC[12](success) -= K4*F2;
1341  fC[13](success) -= K4*F3;
1342  fC[14](success) -= K4*F4;
1343 
1344  return success;
1345 }
TCanvas * c11
TCanvas * c10
TCanvas * c21
static T Abs(const T &x)
Definition: PndCAMath.h:39
TFile * f
Definition: bump_analys.C:12
Double_t z
TCanvas * c20
Double_t y
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 
)

Definition at line 1347 of file PndFTSCATrackParamVector.cxx.

References CAMath::Abs(), c10, c11, c20, c21, FTSCAStripInfo::cos, f, fC, fChi2, fNDF, fP, h1, r, rsqrt(), FTSCAStripInfo::sin, SinPhi(), sqrt(), mrfdata_8b_error::success, and y0.

1350 {
1351  // linearize track in current point, which must be x == x0
1352  // distance between wire and track are : r = - ({x,y,z} - {x0,y0,z0}) * e_o / |e_o|, where ort e_o = |[e_t x e_s]|
1353  const float_v& etx = sqrt( 1 - SinPhi()*SinPhi() );
1354  const float_v& ety = SinPhi();
1355  const float_v& etz = 1.f;
1356  const float_v& eox = - ety * info.cos - etz * info.sin; // e_o[x] = - e_t[y] * e_s[z] + e_t[z] * e_s[y]
1357  const float_v& eoy = etx * info.cos; // e_o[y] = - 0 + e_t[x] * e_s[z]
1358  const float_v& eoz = etx * info.sin; // e_o[z] = - e_t[x] * e_s[y] + 0
1359 
1360  const float_v& iEo = rsqrt( eox*eox + eoy*eoy + eoz*eoz );
1361  const float_v& h0 = eoy*iEo;
1362  const float_v& h1 = eoz*iEo;
1363 
1364  assert( maxSinPhi > 0.f );
1365  float_m success = mask;
1366  success &= ( err2 > 1.e-8f );
1367 
1368  const float_v& c00 = fC[0];
1369  const float_v& c10 = fC[1];
1370  const float_v& c11 = fC[2];
1371 
1372  const float_v zeta = h0*(fP[0] - y0) + h1*(fP[1] - z0) - r;
1373 
1374  // F = CH'
1375  const float_v F0 = h0*c00 + h1*c10;
1376  const float_v F1 = h0*c10 + h1*c11;
1377 
1378  const float_v HCH = ( F0*h0 + F1*h1 );
1379 
1380  const float_v wi = 1/( err2 + HCH );
1381  const float_v zetawi = zeta *wi;
1382 
1383  fChi2(success) += zeta * zetawi;
1384  success &= fChi2 < chi2Cut;
1385  if ( success.isEmpty() ) return success;
1386 
1387  const float_v& c20 = fC[3];
1388  const float_v& c21 = fC[4];
1389  const float_v& c30 = fC[6];
1390  const float_v& c31 = fC[7];
1391  const float_v& c40 = fC[10];
1392  const float_v& c41 = fC[11];
1393 
1394  const float_v F2 = h0*c20 + h1*c21;
1395  const float_v F3 = h0*c30 + h1*c31;
1396  const float_v F4 = h0*c40 + h1*c41;
1397 
1398 
1399  const float_v K1 = F1*wi;
1400  const float_v K2 = F2*wi;
1401  const float_v K3 = F3*wi;
1402  const float_v K4 = F4*wi;
1403 
1404 
1405  const float_v sinPhi = fP[2] - F2*zetawi;
1406 
1407  success &= CAMath::Abs( sinPhi ) < maxSinPhi;
1408 
1409  fNDF( static_cast<int_m>(success) ) += 1; // TODO
1410 
1411  fP[ 0](success) -= F0*zetawi;
1412  fP[ 1](success) -= F1*zetawi;
1413  fP[ 2](success) = sinPhi ;
1414  fP[ 3](success) -= F3*zetawi;
1415  fP[ 4](success) -= F4*zetawi;
1416  fC[ 0](success) -= F0*F0*wi;
1417 
1418  fC[ 1](success) -= K1*F0;
1419  fC[ 2](success) -= K1*F1;
1420 
1421  fC[ 3](success) -= K2*F0;
1422  fC[ 4](success) -= K2*F1;
1423  fC[ 5](success) -= K2*F2;
1424 
1425  fC[ 6](success) -= K3*F0;
1426  fC[ 7](success) -= K3*F1;
1427  fC[ 8](success) -= K3*F2;
1428  fC[ 9](success) -= K3*F3;
1429 
1430  fC[10](success) -= K4*F0;
1431  fC[11](success) -= K4*F1;
1432 
1433  fC[12](success) -= K4*F2;
1434  fC[13](success) -= K4*F3;
1435  fC[14](success) -= K4*F4;
1436 
1437  return success;
1438 }
Double_t z0
Definition: checkhelixhit.C:62
TCanvas * c11
double r
Definition: RiemannTest.C:14
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
TCanvas * c10
TCanvas * c21
static T Abs(const T &x)
Definition: PndCAMath.h:39
Double_t y0
Definition: checkhelixhit.C:71
friend F32vec4 rsqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:32
TFile * f
Definition: bump_analys.C:12
TCanvas * c20
float_v PndFTSCATrackParamVector::GetChi2 ( ) const
inline

Definition at line 1347 of file PndFTSCATrackParamVector.h.

References fChi2.

Referenced by SetTrackParam(), and SetTrackParamOne().

1347 { return fChi2; }
float_v PndFTSCATrackParamVector::GetCosPhi ( ) const
inline

Definition at line 1352 of file PndFTSCATrackParamVector.h.

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

Referenced by GetDCAPoint(), GetS(), Rotate(), and RotateXY().

1352 { return fSignCosPhi*CAMath::Sqrt( float_v( Vc::One ) - SinPhi()*SinPhi() ); }
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
static const fvec One
float_v PndFTSCATrackParamVector::GetCosPhiPositive ( ) const
inline

Definition at line 1351 of file PndFTSCATrackParamVector.h.

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

1351 { return CAMath::Sqrt( float_v( Vc::One ) - SinPhi()*SinPhi() ); }
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
static const fvec One
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

Definition at line 737 of file PndFTSCATrackParamVector.cxx.

References a, CAMath::Abs(), dx, dy, GetCosPhi(), GetDzDs(), GetKappa(), GetS(), GetSinPhi(), GetX(), GetY(), GetZ(), CAMath::Round(), s, sqrt(), CAMath::TwoPi(), x0, and y0.

739 {
740  //* Get the track point closest to the (x,y,z)
741 
742  const float_v &x0 = GetX();
743  const float_v &y0 = GetY();
744  const float_v &k = GetKappa( Bz );
745  const float_v &ex = GetCosPhi();
746  const float_v &ey = GetSinPhi();
747  const float_v &dx = x - x0;
748  const float_v &dy = y - y0;
749  const float_v &ax = dx * k + ey;
750  const float_v &ay = dy * k - ex;
751  const float_v &a = sqrt( ax * ax + ay * ay );
752  *xp = x0 + ( dx - ey * ( ( dx * dx + dy * dy ) * k - 2.f * ( -dx * ey + dy * ex ) ) / ( a + 1.f ) ) / a;
753  *yp = y0 + ( dy + ex * ( ( dx * dx + dy * dy ) * k - 2.f * ( -dx * ey + dy * ex ) ) / ( a + 1.f ) ) / a;
754  const float_v s = GetS( x, y, Bz );
755  *zp = GetZ() + GetDzDs() * s;
756  const float_v dZ = CAMath::Abs( GetDzDs() * CAMath::TwoPi() / k );
757  const float_m mask = CAMath::Abs( k ) > 1.e-2f && dZ > .1f;
758  ( *zp )( mask ) += CAMath::Round( ( z - *zp ) / dZ ) * dZ;
759 }
Double_t x0
Definition: checkhelixhit.C:70
double dy
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
TLorentzVector s
Definition: Pnd2DStar.C:50
static T Round(const T &x)
float_v GetS(const float_v &x, const float_v &y, const float_v &Bz) const
float_v GetKappa(const float_v &Bz) const
static T Abs(const T &x)
Definition: PndCAMath.h:39
Int_t a
Definition: anaLmdDigi.C:126
Double_t y0
Definition: checkhelixhit.C:71
TFile * f
Definition: bump_analys.C:12
Double_t z
double dx
Double_t x
static float TwoPi()
Definition: PndCAMath.h:61
Double_t y
float_v PndFTSCATrackParamVector::GetDist2 ( const PndFTSCATrackParamVector t) const

Definition at line 700 of file PndFTSCATrackParamVector.cxx.

References dx, dy, dz, GetX(), GetY(), and GetZ().

701 {
702  // get squared distance between tracks
703 
704  const float_v &dx = GetX() - t.GetX();
705  const float_v &dy = GetY() - t.GetY();
706  const float_v &dz = GetZ() - t.GetZ();
707  return dx*dx + dy*dy + dz*dz;
708 }
double dy
double dx
float_v PndFTSCATrackParamVector::GetDistXZ2 ( const PndFTSCATrackParamVector t) const

Definition at line 710 of file PndFTSCATrackParamVector.cxx.

References dx, dz, GetX(), and GetZ().

711 {
712  // get squared distance between tracks in X&Z
713 
714  const float_v &dx = GetX() - t.GetX();
715  const float_v &dz = GetZ() - t.GetZ();
716  return dx*dx + dz*dz;
717 }
double dx
float_v PndFTSCATrackParamVector::GetDzDs ( ) const
inline

Definition at line 1344 of file PndFTSCATrackParamVector.h.

References fP.

Referenced by GetDCAPoint().

1344 { return fP[3]; }
float_v PndFTSCATrackParamVector::GetKappa ( const float_v &  Bz) const
inline

Definition at line 1350 of file PndFTSCATrackParamVector.h.

References fP.

Referenced by GetDCAPoint(), and GetS().

1350 { return fP[4]*Bz; }
int_v PndFTSCATrackParamVector::GetNDF ( ) const
inline

Definition at line 1348 of file PndFTSCATrackParamVector.h.

References fNDF.

Referenced by SetTrackParam(), and SetTrackParamOne().

float_v PndFTSCATrackParamVector::GetQPt ( ) const
inline

Definition at line 1345 of file PndFTSCATrackParamVector.h.

References fP.

1345 { return fP[4]; }
float_v PndFTSCATrackParamVector::GetS ( const float_v &  x,
const float_v &  y,
const float_v &  Bz 
) const

Definition at line 720 of file PndFTSCATrackParamVector.cxx.

References CAMath::Abs(), CAMath::ATan2(), f, GetCosPhi(), GetKappa(), GetSinPhi(), GetX(), GetY(), x, and y.

Referenced by GetDCAPoint().

721 {
722  //* Get XY path length to the given point
723 
724  const float_v &k = GetKappa( Bz );
725  const float_v &ex = GetCosPhi();
726  const float_v &ey = GetSinPhi();
727  const float_v &x = _x - GetX();
728  const float_v &y = _y - GetY();
729  float_v dS = x * ex + y * ey;
730  const float_m &mask = CAMath::Abs( k ) > 1.e-4f;
731  if ( !mask.isEmpty() ) {
732  dS( mask ) = CAMath::ATan2( k * dS, 1 + k * ( x * ey - y * ex ) ) / k;
733  }
734  return dS;
735 }
float_v GetKappa(const float_v &Bz) const
static T Abs(const T &x)
Definition: PndCAMath.h:39
static T ATan2(const T &y, const T &x)
TFile * f
Definition: bump_analys.C:12
Double_t x
Double_t y
float_v PndFTSCATrackParamVector::GetSignCosPhi ( ) const
inline

Definition at line 1346 of file PndFTSCATrackParamVector.h.

References fSignCosPhi.

float_v PndFTSCATrackParamVector::GetSinPhi ( ) const
inline

Definition at line 1343 of file PndFTSCATrackParamVector.h.

References fP.

Referenced by GetDCAPoint(), and GetS().

1343 { return fP[2]; }
float_v PndFTSCATrackParamVector::GetX ( ) const
inline

Definition at line 1340 of file PndFTSCATrackParamVector.h.

References fX.

Referenced by GetDCAPoint(), GetDist2(), GetDistXZ2(), and GetS().

1340 { return fX; }
float_v PndFTSCATrackParamVector::GetY ( ) const
inline

Definition at line 1341 of file PndFTSCATrackParamVector.h.

References fP.

Referenced by GetDCAPoint(), GetDist2(), and GetS().

1341 { return fP[0]; }
float_v PndFTSCATrackParamVector::GetZ ( ) const
inline

Definition at line 1342 of file PndFTSCATrackParamVector.h.

References fP.

Referenced by GetDCAPoint(), GetDist2(), and GetDistXZ2().

1342 { return fP[1]; }
void PndFTSCATrackParamVector::InitByHit ( const FTSCAHitV hit,
const PndFTSCAParam param,
const float_v &  dQMom 
)

Definition at line 652 of file PndFTSCATrackParamVector.cxx.

References FTSCAHitV::Angle(), FTSCAHitV::Err2X1(), FTSCAHitV::Err2X2(), SetAngle(), SetChi2(), SetCov(), SetDzDs(), SetNDF(), SetQPt(), SetSignCosPhi(), SetSinPhi(), SetX(), SetY(), SetZ(), FTSCAHitV::X0(), FTSCAHitV::X1(), and FTSCAHitV::X2().

Referenced by PndFTSCAGBTracker::Create1Plets(), and PndFTSCAGBTracker::FitTrackCA().

653 {
654  SetX( hit.X0() );
655  SetY( hit.X1() );
656  SetZ( hit.X2() );
657  SetSinPhi( 0.f );
658  SetSignCosPhi( 1.f );
659  SetDzDs( 0.f );
660  SetQPt( 0.f );
661  SetAngle( hit.Angle() );
662 
663  SetCov( 0, hit.Err2X1() );
664  SetCov( 1, 0.f );
665  SetCov( 2, hit.Err2X2() );
666  SetCov( 3, 0.f );
667  SetCov( 4, 0.f );
668  SetCov( 5, 1.f );
669  SetCov( 6, 0.f );
670  SetCov( 7, 0.f );
671  SetCov( 8, 0.f );
672  SetCov( 9, 1.f );
673  SetCov(10, 0.f );
674  SetCov(11, 0.f );
675  SetCov(12, 0.f );
676  SetCov(13, 0.f );
677  SetCov(14, dQMom*dQMom );
678 
679  SetChi2( 0.f );
680 
681  SetNDF( -3 );
682 
683  UNUSED_PARAM1(param);
684  // const int ista0 = hit.IStation();
685  // const FTSCAStation sta0 = param.Station(ista0);
686 
687  // CAFieldValue B;
688  // sta0.fieldSlice.GetFieldValue( hit.X1(), hit.X2(), B );
689  // SetField( 1, B, sta0.z );
690 
691  // if ( ista0 > 0 ) {
692  // const FTSCAStation staC = param.Station(ista0/2); // The field will be used to propaget to the target. Therefore take station in a middle between target and hit.
693  // const float_v tx1 = hit.X1()/hit.X0(); // suppose target is at 0
694  // const float_v tx2 = hit.X2()/hit.X0(); // suppose target is at 0
695  // staC.fieldSlice.GetFieldValue( tx1*staC.z, tx2*staC.z, B );
696  // SetField( 0, B, staC.z );
697  // }
698 }
float_v X2() const
Definition: FTSCAHitsV.h:51
float_v Err2X1() const
Definition: FTSCAHitsV.h:56
void SetCov(int i, const float_v &v)
void SetSignCosPhi(const float_v &v)
void SetSinPhi(const float_v &v)
float_v X1() const
Definition: FTSCAHitsV.h:50
float_v Angle() const
Definition: FTSCAHitsV.h:73
TFile * f
Definition: bump_analys.C:12
float_v X0() const
Definition: FTSCAHitsV.h:49
float_v Err2X2() const
Definition: FTSCAHitsV.h:58
void PndFTSCATrackParamVector::InitByTarget ( const FTSCATarget target)

Definition at line 619 of file PndFTSCATrackParamVector.cxx.

References FTSCATarget::Err2QMom(), FTSCATarget::Err2X1(), FTSCATarget::Err2X2(), FTSCATarget::NDF(), SetChi2(), SetCov(), SetDzDs(), SetNDF(), SetQPt(), SetSignCosPhi(), SetSinPhi(), SetX(), SetY(), SetZ(), FTSCATarget::X0(), FTSCATarget::X1(), and FTSCATarget::X2().

Referenced by PndFTSCAGBTracker::Create1Plets(), CAFunctionality::FitIteration(), PndFTSCAGBTracker::FitTrackCA(), and PndFTSCAGBTracker::Refit().

620 {
621  SetX( t.X0() );
622  SetY( t.X1() );
623  SetZ( t.X2() );
624  SetSinPhi( 1.f );
625  SetSignCosPhi( 0.f );
626  SetDzDs( 0.f );
627  SetQPt( 0.f );
628  SetCov( 0.f, t.Err2X1() );
629  SetCov( 1, 0.f );
630  SetCov( 2, t.Err2X2() );
631  SetCov( 3, 0.f );
632  SetCov( 4, 0.f );
633  SetCov( 5, 1.f );
634  SetCov( 6, 0.f );
635  SetCov( 7, 0.f );
636  SetCov( 8, 0.f );
637  SetCov( 9, 1.f );
638  SetCov(10, 0.f );
639  SetCov(11, 0.f );
640  SetCov(12, 0.f );
641  SetCov(13, 0.f );
642  SetCov(14, t.Err2QMom() );
643 
644  SetChi2( 0.f );
645 
646  SetNDF( -5 + t.NDF() );
647 
648  // SetField( 0, t.B(), t.X0() );
649 }
void SetCov(int i, const float_v &v)
void SetSignCosPhi(const float_v &v)
void SetSinPhi(const float_v &v)
TFile * f
Definition: bump_analys.C:12
TTree * t
Definition: bump_analys.C:13
void PndFTSCATrackParamVector::InitCovMatrix ( float_v  d2QMom = 0.f)

Definition at line 596 of file PndFTSCATrackParamVector.cxx.

References SetChi2(), SetCov(), and SetNDF().

597 {
598  SetCov( 0.f, 100.f );
599  SetCov( 1, 0.f );
600  SetCov( 2, 100.f );
601  SetCov( 3, 0.f );
602  SetCov( 4, 0.f );
603  SetCov( 5, 1.f );
604  SetCov( 6, 0.f );
605  SetCov( 7, 0.f );
606  SetCov( 8, 0.f );
607  SetCov( 9, 1.f );
608  SetCov(10, 0.f );
609  SetCov(11, 0.f );
610  SetCov(12, 0.f );
611  SetCov(13, 0.f );
612  SetCov(14, d2QMom );
613 
614  SetChi2( 0.f );
615 
616  SetNDF( -5 );
617 }
void SetCov(int i, const float_v &v)
TFile * f
Definition: bump_analys.C:12
void PndFTSCATrackParamVector::InitDirection ( float_v  r0,
float_v  r1,
float_v  r2 
)
inline

Definition at line 1290 of file PndFTSCATrackParamVector.h.

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

Referenced by PndFTSCAGBTracker::Create1Plets(), CAFunctionality::FitIteration(), PndFTSCAGBTracker::FitTrackCA(), and PndFTSCAGBTracker::Refit().

1291  {
1292  const float_v r = sqrt( r0*r0+r1*r1 );
1293  SetSinPhi( r1/r );
1294  SetSignCosPhi( r0/abs(r0) );
1295  SetDzDs( r2/r );
1296  }
double r
Definition: RiemannTest.C:14
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
void SetSignCosPhi(const float_v &v)
void SetSinPhi(const float_v &v)
double r1
double r2
int_v PndFTSCATrackParamVector::NDF ( ) const
inline
const float_v& PndFTSCATrackParamVector::Par ( int  i) const
inline

Definition at line 1360 of file PndFTSCATrackParamVector.h.

References fP, and i.

Referenced by PndFTSCAGBTracker::FitTrack(), operator<<(), SetTrackParam(), and SetTrackParamOne().

1360 { return fP[i]; }
Int_t i
Definition: run_full.C:25
const float_v* PndFTSCATrackParamVector::Par ( ) const
inlineprivate

Definition at line 1365 of file PndFTSCATrackParamVector.h.

References fP.

Referenced by ConvertTrackParamToVector().

1365 { return fP; }
float_v* PndFTSCATrackParamVector::Par ( )
inlineprivate

Definition at line 1367 of file PndFTSCATrackParamVector.h.

References fP.

1367 { return fP; }
float_v PndFTSCATrackParamVector::QMomentum ( ) const
inline

Definition at line 1324 of file PndFTSCATrackParamVector.h.

References QPt().

Referenced by CAFunctionality::Fit(), and CAFunctionality::FitIteration().

1324 { return QPt(); } // used for triplets comparison
float_v PndFTSCATrackParamVector::QPt ( ) const
inline
float_m PndFTSCATrackParamVector::Rotate ( const float_v &  alpha,
PndFTSCATrackLinearisationVector t0,
const float  maxSinPhi = .999f,
const float_m &  mask = float_m( true ) 
)

Definition at line 1116 of file PndFTSCATrackParamVector.cxx.

References CAMath::Abs(), alpha, CAMath::Cos(), PndFTSCATrackLinearisationVector::CosPhi(), d, f, fAlpha, fC, PndFTSCATrackLinearisationVector::SetCosPhi(), PndFTSCATrackLinearisationVector::SetSinPhi(), SetSinPhi(), SetX(), SetY(), CAMath::Sin(), PndFTSCATrackLinearisationVector::SinPhi(), SinPhi(), X(), x0, Y(), and y0.

Referenced by PndFTSCAGBTracker::FitTrack(), and Transport().

1118 {
1119  //* Rotate the coordinate system in XY on the angle alpha
1120  if ( (CAMath::Abs(alpha) < 1e-6f || !mask).isFull() ) return mask;
1121 
1122  const float_v cA = CAMath::Cos( alpha );
1123  const float_v sA = CAMath::Sin( alpha );
1124  const float_v x0 = X(), y0 = Y(), sP = t0.SinPhi(), cP = t0.CosPhi();
1125  const float_v cosPhi = cP * cA + sP * sA;
1126  const float_v sinPhi = -cP * sA + sP * cA;
1127 
1128  float_m ReturnMask(mask);
1129  ReturnMask &= (!( CAMath::Abs( sinPhi ) > maxSinPhi || CAMath::Abs( cosPhi ) < 1.e-4f || CAMath::Abs( cP ) < 1.e-4f ));
1130 
1131  float_v tmp = alpha*0.15915494f;// 1/(2.f*3.1415f);
1132  ReturnMask &= abs(tmp - round(tmp)) < 0.167f; // allow turn by 60 degree only TODO scalar
1133 
1134  //float J[5][5] = { { j0, 0, 0, 0, 0 }, // Y
1135  // { 0, 1, 0, 0, 0 }, // Z
1136  // { 0, 0, j2, 0, 0 }, // SinPhi
1137  // { 0, 0, 0, 1, 0 }, // DzDs
1138  // { 0, 0, 0, 0, 1 } }; // Kappa
1139 
1140  const float_v j0 = cP / cosPhi;
1141  const float_v j2 = cosPhi / cP;
1142  const float_v d = SinPhi() - sP;
1143 
1144  SetX( x0*cA + y0*sA, ReturnMask );
1145  SetY(-x0*sA + y0*cA, ReturnMask );
1146  t0.SetCosPhi( cosPhi );
1147  t0.SetSinPhi( sinPhi );
1148 
1149  SetSinPhi( sinPhi + j2*d, ReturnMask );
1150 
1151  fC[0](ReturnMask) *= j0 * j0;
1152  fC[1](ReturnMask) *= j0;
1153  fC[3](ReturnMask) *= j0;
1154  fC[6](ReturnMask) *= j0;
1155  fC[10](ReturnMask) *= j0;
1156 
1157  fC[3](ReturnMask) *= j2;
1158  fC[4](ReturnMask) *= j2;
1159  fC[5](ReturnMask) *= j2 * j2;
1160  fC[8](ReturnMask) *= j2;
1161  fC[12](ReturnMask) *= j2;
1162 
1163  fAlpha(ReturnMask) += alpha;
1164 
1165  return ReturnMask;
1166 }
Double_t x0
Definition: checkhelixhit.C:70
TObjArray * d
static T Sin(const T &x)
Definition: PndCAMath.h:42
void SetSinPhi(const float_v &v)
static T Cos(const T &x)
Definition: PndCAMath.h:43
static T Abs(const T &x)
Definition: PndCAMath.h:39
Double_t y0
Definition: checkhelixhit.C:71
TFile * f
Definition: bump_analys.C:12
double alpha
Definition: f_Init.h:9
float_m PndFTSCATrackParamVector::Rotate ( const float_v &  alpha,
const float  maxSinPhi = .999f,
const float_m &  mask = float_m( true ) 
)
inline

Definition at line 1855 of file PndFTSCATrackParamVector.h.

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

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

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

1904 {
1905  //* Rotate the coordinate system in XY on the angle alpha
1906  if ( (abs(alpha) < 1e-6f || !mask).isFull() ) return;
1907 
1908  const float_v cA = CAMath::Cos( alpha );
1909  const float_v sA = CAMath::Sin( alpha );
1910 
1911  x(mask) = ( X()*cA + Y()*sA );
1912  y(mask) = ( -X()*sA + Y()*cA );
1913  sin(mask) = -GetCosPhi() * sA + SinPhi() * cA;
1914 }
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
static T Sin(const T &x)
Definition: PndCAMath.h:42
static T Cos(const T &x)
Definition: PndCAMath.h:43
TFile * f
Definition: bump_analys.C:12
Double_t x
Double_t y
double alpha
Definition: f_Init.h:9
void PndFTSCATrackParamVector::SetAngle ( const float_v &  v)
inline
void PndFTSCATrackParamVector::SetAngle ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1420 of file PndFTSCATrackParamVector.h.

References fAlpha, and v.

1420 { fAlpha(m) = v; }
__m128 m
Definition: P4_F32vec4.h:28
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetChi2 ( const float_v &  v)
inline
void PndFTSCATrackParamVector::SetChi2 ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1414 of file PndFTSCATrackParamVector.h.

References fChi2, and v.

1414 { fChi2(m) = v; }
__m128 m
Definition: P4_F32vec4.h:28
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetCov ( int  i,
const float_v &  v 
)
inline
void PndFTSCATrackParamVector::SetCov ( int  i,
const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1396 of file PndFTSCATrackParamVector.h.

References i, and m.

1396 { fC[i]( m ) = v; }
Int_t i
Definition: run_full.C:25
__m128 m
Definition: P4_F32vec4.h:28
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetDzDs ( const float_v &  v)
inline

Definition at line 1406 of file PndFTSCATrackParamVector.h.

References fP, and v.

Referenced by InitByHit(), InitByTarget(), and InitDirection().

1406 { fP[3] = v; }
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetDzDs ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1407 of file PndFTSCATrackParamVector.h.

References fP, and m.

1407 { fP[3]( m ) = v; }
__m128 m
Definition: P4_F32vec4.h:28
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetErr2QPt ( float_v  v)
inline

Definition at line 1424 of file PndFTSCATrackParamVector.h.

References v.

1424 { fC[14] = v; }
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetErr2Y ( float_v  v)
inline

Definition at line 1422 of file PndFTSCATrackParamVector.h.

References v.

1422 { fC[0] = v; }
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetErr2Z ( float_v  v)
inline

Definition at line 1423 of file PndFTSCATrackParamVector.h.

References v.

1423 { fC[2] = v; }
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetNDF ( int  v)
inline
void PndFTSCATrackParamVector::SetNDF ( const int_v &  v)
inline

Definition at line 1416 of file PndFTSCATrackParamVector.h.

References fNDF, and v.

1416 { fNDF = v; }
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetNDF ( const int_v &  v,
const int_m &  m 
)
inline

Definition at line 1417 of file PndFTSCATrackParamVector.h.

References fNDF, and v.

1417 { fNDF(m) = v; }
__m128 m
Definition: P4_F32vec4.h:28
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetPar ( int  i,
const float_v &  v 
)
inline

Definition at line 1393 of file PndFTSCATrackParamVector.h.

References fP, i, and v.

Referenced by ConvertTrackParamToVector().

1393 { fP[i] = v; }
Int_t i
Definition: run_full.C:25
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetPar ( int  i,
const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1394 of file PndFTSCATrackParamVector.h.

References fP, i, and m.

1394 { fP[i]( m ) = v; }
Int_t i
Definition: run_full.C:25
__m128 m
Definition: P4_F32vec4.h:28
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetQMomentum ( const float_v &  v)
inline

Definition at line 1409 of file PndFTSCATrackParamVector.h.

References SetQPt().

Referenced by CAFunctionality::FitIteration().

1409 { SetQPt(v); }
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetQPt ( const float_v &  v)
inline

Definition at line 1408 of file PndFTSCATrackParamVector.h.

References fP, and v.

Referenced by PndFTSCAGBTracker::FitTrack(), PndFTSCAGBTracker::FitTrackCA(), InitByHit(), InitByTarget(), and SetQMomentum().

1408 { fP[4] = v; }
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetQPt ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1410 of file PndFTSCATrackParamVector.h.

References fP, and m.

1410 { fP[4]( m ) = v; }
__m128 m
Definition: P4_F32vec4.h:28
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetSignCosPhi ( const float_v &  v)
inline
void PndFTSCATrackParamVector::SetSignCosPhi ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1412 of file PndFTSCATrackParamVector.h.

References fSignCosPhi, and v.

1412 { fSignCosPhi(m) = v; }
__m128 m
Definition: P4_F32vec4.h:28
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetSinPhi ( const float_v &  v)
inline

Definition at line 1404 of file PndFTSCATrackParamVector.h.

References fP, and v.

Referenced by InitByHit(), InitByTarget(), InitDirection(), and Rotate().

1404 { fP[2] = v; }
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetSinPhi ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1405 of file PndFTSCATrackParamVector.h.

References fP, and m.

1405 { fP[2]( m ) = v; }
__m128 m
Definition: P4_F32vec4.h:28
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetTrackParam ( const PndFTSCATrackParamVector param,
const float_m &  m = float_m( true ) 
)
inline

Definition at line 1371 of file PndFTSCATrackParamVector.h.

References Angle(), Cov(), fAlpha, fChi2, fNDF, fP, fSignCosPhi, fX, GetChi2(), GetNDF(), i, m, Par(), SignCosPhi(), and X().

Referenced by PndFTSCAGBTracker::FitTrackCA().

1372  {
1373  for(int i=0; i<5; i++) fP[i](m) = param.Par()[i];
1374  for(int i=0; i<15; i++) fC[i](m) = param.Cov()[i];
1375  fX(m) = param.X();
1376  fSignCosPhi(m) = param.SignCosPhi();
1377  fChi2(m) = param.GetChi2();
1378  fNDF(static_cast<int_m>(m)) = param.GetNDF();
1379  fAlpha(static_cast<int_m>(m)) = param.Angle();
1380  }
Int_t i
Definition: run_full.C:25
__m128 m
Definition: P4_F32vec4.h:28
const float_v & Cov(int i) const
const float_v & Par(int i) const
void PndFTSCATrackParamVector::SetTrackParamOne ( int  iV,
const PndFTSCATrackParamVector param,
int  iVa 
)
inline

Definition at line 1382 of file PndFTSCATrackParamVector.h.

References Angle(), Cov(), fAlpha, fChi2, fNDF, fP, fSignCosPhi, fX, GetChi2(), GetNDF(), i, Par(), SignCosPhi(), and X().

Referenced by FTSCANPletV::CopyOne().

1383  {
1384  for(int i=0; i<5; i++) fP[i][iV] = param.Par()[i][iVa];
1385  for(int i=0; i<15; i++) fC[i][iV] = param.Cov()[i][iVa];
1386  fX[iV] = param.X()[iVa];
1387  fSignCosPhi[iV] = param.SignCosPhi()[iVa];
1388  fChi2[iV] = param.GetChi2()[iVa];
1389  fNDF[iV] = param.GetNDF()[iVa];
1390  fAlpha[iV] = param.Angle()[iVa];
1391  }
Int_t i
Definition: run_full.C:25
const float_v & Cov(int i) const
const float_v & Par(int i) const
void PndFTSCATrackParamVector::SetX ( const float_v &  v)
inline
void PndFTSCATrackParamVector::SetX ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1401 of file PndFTSCATrackParamVector.h.

References fX, and v.

1401 { fX( m ) = v; }
__m128 m
Definition: P4_F32vec4.h:28
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetY ( const float_v &  v)
inline

Definition at line 1399 of file PndFTSCATrackParamVector.h.

References fP, and v.

Referenced by PndFTSCAGBTracker::FitTrackCA(), InitByHit(), InitByTarget(), PndFTSCAGBTracker::Refit_1(), and Rotate().

1399 { fP[0] = v; }
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetY ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1402 of file PndFTSCATrackParamVector.h.

References fP, and m.

1402 { fP[0]( m ) = v; }
__m128 m
Definition: P4_F32vec4.h:28
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetZ ( const float_v &  v)
inline

Definition at line 1400 of file PndFTSCATrackParamVector.h.

References fP, and v.

Referenced by PndFTSCAGBTracker::FitTrackCA(), InitByHit(), and InitByTarget().

1400 { fP[1] = v; }
__m128 v
Definition: P4_F32vec4.h:4
void PndFTSCATrackParamVector::SetZ ( const float_v &  v,
const float_m &  m 
)
inline

Definition at line 1403 of file PndFTSCATrackParamVector.h.

References fP, and m.

1403 { fP[1]( m ) = v; }
__m128 m
Definition: P4_F32vec4.h:28
__m128 v
Definition: P4_F32vec4.h:4
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 1330 of file PndFTSCATrackParamVector.h.

References fSignCosPhi.

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

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

Definition at line 1440 of file PndFTSCATrackParamVector.cxx.

References CalculateFitParameters(), PndFTSCAParam::cBz(), PndFTSCAParam::GetX0(), PndFTSCAParam::GetXOverX0(), PndFTSCAParam::GetXTimesRho(), and TransportToX0WithMaterial().

Referenced by PndFTSCAGBTracker::Create1Plets(), CAFunctionality::FitIteration(), PndFTSCAGBTracker::FitTrackCA(), PndFTSCAGBTracker::PickUpHits(), and PndFTSCAGBTracker::Refit().

1441 {
1442  float_m active = mask;
1443 
1444  PndFTSCATrackFitParam fitPar;
1445  CalculateFitParameters( fitPar );
1447  active &= TransportToX0WithMaterial( param.GetX0( ista, active ), tE, fitPar, param.GetXOverX0(int_v(ista),active), param.GetXTimesRho(int_v(ista),active), param.cBz(), 0.999f, active );
1448  return active;
1449 }
float GetX0(short iSt) const
Definition: PndFTSCAParam.h:61
float GetXOverX0(short iSt) const
Definition: PndFTSCAParam.h:65
float cBz() const
Definition: PndFTSCAParam.h:49
void CalculateFitParameters(PndFTSCATrackFitParam &par, const float_v &mass=0.13957f)
float_m TransportToX0WithMaterial(const float_v &x, const float_v &XOverX0, const float_v &XThimesRho, const float_v &Bz, const float maxSinPhi=.999f)
float GetXTimesRho(short iSt) const
Definition: PndFTSCAParam.h:66
float_m PndFTSCATrackParamVector::Transport ( const FTSCAHitV hit,
const PndFTSCAParam p,
const float_m &  mask = float_m( true ) 
)

Definition at line 1451 of file PndFTSCATrackParamVector.cxx.

References CAMath::Abs(), FTSCAHitV::Angle(), CalculateFitParameters(), PndFTSCAParam::cBz(), fAlpha, fX, PndFTSCAParam::GetXOverX0(), PndFTSCAParam::GetXTimesRho(), FTSCAHitV::IStations(), FTSCAHitV::IsValid(), Rotate(), TransportToX0WithMaterial(), and FTSCAHitV::X0().

1452 {
1453  if ( ((CAMath::Abs(-fAlpha + hit.Angle()) < 1e-6f && CAMath::Abs(-fX + hit.X0()) < 2e-4f) || !mask).isFull() ) return mask;
1454  float_m active = mask & hit.IsValid();
1456  active &= Rotate( -fAlpha + hit.Angle(), tR, .999f, active );
1457 
1458  PndFTSCATrackFitParam fitPar;
1459  CalculateFitParameters( fitPar );
1460  active &= TransportToX0WithMaterial( hit.X0(), tR, fitPar, param.GetXOverX0(hit.IStations(),active), param.GetXTimesRho(hit.IStations(),active), param.cBz(), 0.999f, active );
1461  return active;
1462 }
float_m IsValid() const
Definition: FTSCAHitsV.h:38
float_v Angle() const
Definition: FTSCAHitsV.h:73
float_m Rotate(const float_v &alpha, PndFTSCATrackLinearisationVector &t0, const float maxSinPhi=.999f, const float_m &mask=float_m(true))
int_v IStations() const
Definition: FTSCAHitsV.h:45
static T Abs(const T &x)
Definition: PndCAMath.h:39
TFile * f
Definition: bump_analys.C:12
void CalculateFitParameters(PndFTSCATrackFitParam &par, const float_v &mass=0.13957f)
float_m TransportToX0WithMaterial(const float_v &x, const float_v &XOverX0, const float_v &XThimesRho, const float_v &Bz, const float maxSinPhi=.999f)
float_v X0() const
Definition: FTSCAHitsV.h:49
float_m PndFTSCATrackParamVector::Transport ( const FTSCAHit hit,
const PndFTSCAParam p,
const float_m &  mask = float_m( true ) 
)

Definition at line 1477 of file PndFTSCATrackParamVector.cxx.

References CAMath::Abs(), FTSCAHit::Angle(), CalculateFitParameters(), PndFTSCAParam::cBz(), fAlpha, fX, PndFTSCAParam::GetXOverX0(), PndFTSCAParam::GetXTimesRho(), FTSCAHit::IStation(), Rotate(), TransportToX0WithMaterial(), and FTSCAHit::X0().

1478 {
1479  if ( ((CAMath::Abs(-fAlpha + hit.Angle()) < 1e-6f && CAMath::Abs(-fX + hit.X0()) < 2e-4f) || !mask).isFull() ) return mask;
1480 
1481  float_m active = mask;
1482 
1484  active &= Rotate( -fAlpha + hit.Angle(), tR, .999f, active );
1485 
1486  PndFTSCATrackFitParam fitPar;
1487  CalculateFitParameters( fitPar );
1488  active &= TransportToX0WithMaterial( hit.X0(), tR, fitPar, param.GetXOverX0(hit.IStation()), param.GetXTimesRho(hit.IStation()), param.cBz(), 0.999f, active );
1489  return active;
1490 }
char IStation() const
Definition: FTSCAHits.h:33
float Angle() const
Definition: FTSCAHits.h:71
float_m Rotate(const float_v &alpha, PndFTSCATrackLinearisationVector &t0, const float maxSinPhi=.999f, const float_m &mask=float_m(true))
float X0() const
Definition: FTSCAHits.h:41
static T Abs(const T &x)
Definition: PndCAMath.h:39
TFile * f
Definition: bump_analys.C:12
void CalculateFitParameters(PndFTSCATrackFitParam &par, const float_v &mass=0.13957f)
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 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 1994 of file PndFTSCATrackParamVector.h.

References CAMath::Abs(), dx, DzDs(), h2, h4, One, CAMath::RSqrt(), SinPhi(), and X().

Referenced by AddTarget().

1995 {
1996  const float_v &ey = SinPhi();
1997  const float_v &dx = x - X();
1998  const float_v &exi = CAMath::RSqrt( float_v( Vc::One ) - ey * ey ); // RSqrt
1999 
2000  const float_v &dxBz = dx * cBz;
2001  const float_v &dS = dx * exi;
2002  const float_v &h2 = dS * exi * exi;
2003  const float_v &h4 = .5f * h2 * dxBz;
2004 
2005  float_m mask = active && CAMath::Abs( exi ) <= 1.e4f;
2006 
2007  extrDy(active) = dS*ey;
2008  extrDz(active) = dS*DzDs();
2009  J[0](active) = dS;
2010  J[1](active) = 0;
2011  J[2](active) = h4;
2012  J[3](active) = dS;
2013  J[4](active) = dS;
2014  J[5](active) = 0;
2015  return active;
2016 }
TH1F * h4
static T Abs(const T &x)
Definition: PndCAMath.h:39
static T RSqrt(const T &x)
Definition: PndCAMath.h:38
double dx
Double_t x
static const fvec One
float_m PndFTSCATrackParamVector::TransportToX0 ( const float_v &  x,
const float_v &  Bz,
const float  maxSinPhi = .999f,
const float_m &  mask = float_m( true ) 
)

Definition at line 880 of file PndFTSCATrackParamVector.cxx.

References t0.

Referenced by TransportToX0WithMaterial().

882 {
883  //* Transport the track parameters to X=x
884 
885 // assert( ( x == 0 && mask ).isEmpty() );
886 
888 
889  return TransportToX0( x, t0, Bz, maxSinPhi, 0, mask );
890 }
Int_t t0
Definition: hist-t7.C:106
Double_t x
float_m TransportToX0(const float_v &x, const float_v &Bz, const float maxSinPhi=.999f, const float_m &mask=float_m(true))
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 ) 
)

Definition at line 767 of file PndFTSCATrackParamVector.cxx.

References CAMath::Abs(), CAMath::ASin(), c10, c11, c20, c21, c22, PndFTSCATrackLinearisationVector::CosPhi(), d, debugKF(), dx, dy, PndFTSCATrackLinearisationVector::DzDs(), f, fC, fP, fX, h2, h4, CAMath::Max(), CAMath::Min(), One, PndFTSCATrackLinearisationVector::QPt(), PndFTSCATrackLinearisationVector::SetCosPhi(), PndFTSCATrackLinearisationVector::SetSinPhi(), PndFTSCATrackLinearisationVector::SinPhi(), CAMath::Sqrt(), X(), Y(), Z(), and Zero.

768 {
769  //* Transport the track parameters to X=x, using linearization at t0, and the field value Bz
770  //* maxSinPhi is the max. allowed value for |t0.SinPhi()|
771  //* linearisation of trajectory t0 is also transported to X=x,
772  //* returns 1 if OK
773  //*
774 // std::cout << t0.QPt()[0] << " " << fP[4][0] << std::endl;
775  const float_v ex = t0.CosPhi();
776  const float_v ey = t0.SinPhi();
777  const float_v k = t0.QPt() * Bz;
778  const float_v dx = x - X();
779 
780  float_v ey1 = k * dx + ey;
781 
782  // check for intersection with X=x
783 
784  float_v ex1 = CAMath::Sqrt( 1.f - ey1 * ey1 );
785  ex1( ex < Vc::Zero ) = -ex1;
786 
787  const float_v dx2 = dx * dx;
788  const float_v ss = ey + ey1;
789  const float_v cc = ex + ex1;
790 
791  const float_v cci = 1.f / cc;
792  const float_v exi = 1.f / ex;
793  const float_v ex1i = 1.f / ex1;
794 
795  float_m mask = _mask && CAMath::Abs( ey1 ) <= maxSinPhi && CAMath::Abs( cc ) >= 1.e-4f && CAMath::Abs( ex ) >= 1.e-4f && CAMath::Abs( ex1 ) >= 1.e-4f;
796 
797  const float_v tg = ss * cci; // tan((phi1+phi)/2)
798 
799  const float_v dy = dx * tg;
800  float_v dl = dx * CAMath::Sqrt( 1.f + tg * tg );
801 
802  dl( cc < Vc::Zero ) = -dl;
803  const float_v dSin = CAMath::Max( float_v( -1.f ), CAMath::Min( float_v( Vc::One ), dl * k * 0.5f ) );
804  const float_v dS = ( CAMath::Abs( k ) > 1.e-4f ) ? ( 2 * CAMath::ASin( dSin ) / k ) : dl;
805  const float_v dz = dS * t0.DzDs();
806 
807  if ( DL ) {
808  ( *DL )( mask ) = -dS * CAMath::Sqrt( 1.f + t0.DzDs() * t0.DzDs() );
809  }
810 
811  const float_v d[3] = { fP[2] - t0.SinPhi(), fP[3] - t0.DzDs(), fP[4] - t0.QPt() };
812 
813  //float H0[5] = { 1,0, h2, 0, h4 };
814  //float H1[5] = { 0, 1, 0, dS, 0 };
815  //float H2[5] = { 0, 0, 1, 0, dxBz };
816  //float H3[5] = { 0, 0, 0, 1, 0 };
817  //float H4[5] = { 0, 0, 0, 0, 1 };
818 
819  const float_v h2 = dx * ( 1.f + ey * ey1 + ex * ex1 ) * exi * ex1i * cci;
820  const float_v h4 = dx2 * ( cc + ss * ey1 * ex1i ) * cci * cci * Bz;
821  const float_v dxBz = dx * Bz;
822 
823  ex1( !mask ) = ex;
824  ey1( !mask ) = ey;
825 
826 // std::cout << "QPt " << t0.QPt()[0] << " " << fP[4][0] <<
827 // " CosPhi " << t0.CosPhi()[0] << " " << ex1[0] << std::endl;
828 
829  t0.SetCosPhi( ex1 );
830  t0.SetSinPhi( ey1 );
831 
832  fX( mask ) = X() + dx;
833  fP[0]( mask ) = Y() + dy + h2 * d[0] + h4 * d[2];
834  fP[1]( mask ) = Z() + dz + dS * d[1];
835  fP[2]( mask ) = t0.SinPhi() + d[0] + dxBz * d[2];
836  fP[2]( CAMath::Abs(fP[2]) > maxSinPhi ) = t0.SinPhi();
837 
838  const float_v c00 = fC[0];
839  const float_v c10 = fC[1];
840  const float_v c11 = fC[2];
841  const float_v c20 = fC[3];
842  const float_v c21 = fC[4];
843  const float_v c22 = fC[5];
844  const float_v c30 = fC[6];
845  const float_v c31 = fC[7];
846  const float_v c32 = fC[8];
847  const float_v c33 = fC[9];
848  const float_v c40 = fC[10];
849  const float_v c41 = fC[11];
850  const float_v c42 = fC[12];
851  const float_v c43 = fC[13];
852  const float_v c44 = fC[14];
853 
854  fC[0]( mask ) = ( c00 + h2 * h2 * c22 + h4 * h4 * c44
855  + 2.f * ( h2 * c20 + h4 * c40 + h2 * h4 * c42 ) );
856 
857  fC[1]( mask ) = c10 + h2 * c21 + h4 * c41 + dS * ( c30 + h2 * c32 + h4 * c43 );
858  fC[2]( mask ) = c11 + 2.f * dS * c31 + dS * dS * c33;
859 
860  fC[3]( mask ) = c20 + h2 * c22 + h4 * c42 + dxBz * ( c40 + h2 * c42 + h4 * c44 );
861  fC[4]( mask ) = c21 + dS * c32 + dxBz * ( c41 + dS * c43 );
862  fC[5]( mask ) = c22 + 2.f * dxBz * c42 + dxBz * dxBz * c44;
863 
864  fC[6]( mask ) = c30 + h2 * c32 + h4 * c43;
865  fC[7]( mask ) = c31 + dS * c33;
866  fC[8]( mask ) = c32 + dxBz * c43;
867  fC[9]( mask ) = c33;
868 
869  fC[10]( mask ) = c40 + h2 * c42 + h4 * c44;
870  fC[11]( mask ) = c41 + dS * c43;
871  fC[12]( mask ) = c42 + dxBz * c44;
872  fC[13]( mask ) = c43;
873  fC[14]( mask ) = c44;
874 //std::cout << fC[10] << " " << fC[11] <<" " << h2 << " " << h4 << " " << c44 << " " << mask << " " << Bz << std::endl;
875  debugKF() << mask << "\n" << *this << std::endl;
876  return mask;
877 }
static T ASin(const T &x)
double dy
TCanvas * c11
TObjArray * d
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
TH1F * h4
TCanvas * c10
TCanvas * c21
static const fvec Zero
static T Abs(const T &x)
Definition: PndCAMath.h:39
static T Min(const T &x, const T &y)
Definition: PndCAMath.h:35
TFile * f
Definition: bump_analys.C:12
double dx
TCanvas * c22
TCanvas * c20
static T Max(const T &x, const T &y)
Definition: PndCAMath.h:36
Double_t x
static const fvec One
PndFTSCANoDebugStream & debugKF()
Definition: debug.h:74
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(), c20, c22, debugKF(), dx, DzDs(), fP, fX, h2, h4, One, QPt(), CAMath::RSqrt(), SinPhi(), and X().

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

Definition at line 929 of file PndFTSCATrackParamVector.cxx.

References CalculateFitParameters(), and par.

Referenced by PndFTSCAGBTracker::FitTrack(), Transport(), and TransportToX0WithMaterial().

931 {
932  //* Transport the track parameters to X=x taking into account material budget
933 
934  PndFTSCATrackFitParam par;
935  CalculateFitParameters( par );
936  return TransportToX0WithMaterial( x, par, XOverX0, XThimesRho, Bz, maxSinPhi );
937 }
Double_t par[3]
void CalculateFitParameters(PndFTSCATrackFitParam &par, const float_v &mass=0.13957f)
float_m TransportToX0WithMaterial(const float_v &x, const float_v &XOverX0, const float_v &XThimesRho, const float_v &Bz, const float maxSinPhi=.999f)
Double_t x
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 ) 
)

Definition at line 894 of file PndFTSCATrackParamVector.cxx.

References CorrectForMeanMaterial(), PndFTSCATrackLinearisationVector::DzDs(), PndFTSCATrackLinearisationVector::SinPhi(), sqrt(), and TransportToX0().

896 {
897  //* Transport the track parameters to X=x taking into account material budget
898 // UNUSED_PARAM1(par);
899  // const float kRho = 1.025e-3f;//0.9e-3;
900 // const float kRadLen = 29.532f;//28.94;
901 // const float kRhoOverRadLen = kRho / kRadLen;
902  // const float kRhoOverRadLen = 7.68e-5;
903  float_v dl;
904 
905  float_m mask = mask_ && TransportToX0( x, t0, Bz, maxSinPhi, &dl, mask_ );
906 // if ( !mask.isEmpty() ) {
907 // CorrectForMeanMaterial( dl * kRhoOverRadLen, dl * kRho, par, mask );
908 // }
909 
910  const float_v dzds2 = t0.DzDs()*t0.DzDs();
911  const float_v sphi2 = t0.SinPhi()*t0.SinPhi();
912  float_v koeff = sqrt( (1 + dzds2)*( 1 + sphi2 ) ); // = dr/dx
913 // mask &= CorrectForMeanMaterial( XOverX0*koeff, XThimesRho*koeff, par, mask );
914  CorrectForMeanMaterial( XOverX0*koeff, XThimesRho*koeff, par, mask );
915 
916  return mask;
917 }
float_m CorrectForMeanMaterial(const float_v &xOverX0, const float_v &xTimesRho, const PndFTSCATrackFitParam &par, const float_m &_mask)
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
Double_t par[3]
Double_t x
float_m TransportToX0(const float_v &x, const float_v &Bz, const float maxSinPhi=.999f, const float_m &mask=float_m(true))
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 
)

Definition at line 920 of file PndFTSCATrackParamVector.cxx.

References t0, and TransportToX0WithMaterial().

922 {
923  //* Transport the track parameters to X=x taking into account material budget
924 
926  return TransportToX0WithMaterial( x, t0, par, XOverX0, XThimesRho, Bz, maxSinPhi );
927 }
Double_t par[3]
Int_t t0
Definition: hist-t7.C:106
float_m TransportToX0WithMaterial(const float_v &x, const float_v &XOverX0, const float_v &XThimesRho, const float_v &Bz, const float maxSinPhi=.999f)
Double_t x
float_v PndFTSCATrackParamVector::Tx1 ( ) const
inline

Definition at line 1322 of file PndFTSCATrackParamVector.h.

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

1322 { return SinPhi()/(SignCosPhi()*sqrt( 1 - SinPhi()*SinPhi() )); } // CHECKME
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
float_v PndFTSCATrackParamVector::Tx2 ( ) const
inline

Definition at line 1323 of file PndFTSCATrackParamVector.h.

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

1323 { return DzDs()/(SignCosPhi()*sqrt( 1 - SinPhi()*SinPhi() )); } // dx2/dx0 = dz/dx
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
float_v PndFTSCATrackParamVector::X ( ) const
inline
float_v PndFTSCATrackParamVector::X0 ( ) const
inline

Definition at line 1319 of file PndFTSCATrackParamVector.h.

References X().

Referenced by AddTarget().

1319 { return X(); }
float_v PndFTSCATrackParamVector::X1 ( ) const
inline

Definition at line 1320 of file PndFTSCATrackParamVector.h.

References Y().

Referenced by PndFTSCAGBTracker::PickUpHits().

1320 { return Y(); }
float_v PndFTSCATrackParamVector::X2 ( ) const
inline

Definition at line 1321 of file PndFTSCATrackParamVector.h.

References Z().

1321 { return Z(); }
float_v PndFTSCATrackParamVector::Y ( ) const
inline
float_v PndFTSCATrackParamVector::Z ( ) const
inline

Friends And Related Function Documentation

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

Definition at line 1543 of file PndFTSCATrackParamVector.cxx.

1544 {
1545  if ( out == std::cerr ) {
1546  out << "------------------------------ Track Param ------------------------------"
1547  << "\n X: " << t.X()
1548  << "\n SignCosPhi: " << t.SignCosPhi()
1549  << "\n Chi2: " << t.Chi2()
1550  << "\n NDF: " << t.NDF()
1551  << "\n Y: " << t.Par()[0]
1552  << "\n Z: " << t.Par()[1]
1553  << "\n SinPhi: " << t.Par()[2]
1554  << "\n DzDs: " << t.Par()[3]
1555  << "\n q/Pt: " << t.Par()[4]
1556  << "\nCovariance Matrix\n";
1557  int i = 0;
1558  out << std::setprecision( 2 );
1559  for ( int step = 1; step <= 5; ++step ) {
1560  int end = i + step;
1561  for ( ; i < end; ++i ) {
1562  out << t.Cov()[i] << '\t';
1563  }
1564  out << "\n";
1565  }
1566  out << std::setprecision( 6 );
1567  return out << std::endl;
1568  }
1569  for ( int j = 0; j < uint_v::Size; ++j ) {
1570  out << t.X()[j] << " "
1571  << t.SignCosPhi()[j] << " "
1572  << t.Chi2()[j] << " "
1573  << t.NDF()[j]
1574  << std::endl;
1575  for ( int i = 0; i < 5; i++ ) out << t.Par()[i][j] << " ";
1576  out << std::endl;
1577  for ( int i = 0; i < 15; i++ ) out << t.Cov()[i][j] << " ";
1578  out << std::endl;
1579  }
1580  return out;
1581 }
Int_t i
Definition: run_full.C:25
TFile * out
Definition: reco_muo.C:20
const float_v & Cov(int i) const
const float_v & Par(int i) const
std::istream& operator>> ( std::istream ,
PndFTSCATrackParamVector  
)
friend

Definition at line 1522 of file PndFTSCATrackParamVector.cxx.

1523 {
1524  float_v::Memory x, s, p[5], c[15], chi2;
1525  int_v::Memory ndf;
1526  for ( int j = 0; j < uint_v::Size; ++j ) {
1527  in >> x[j];
1528  in >> s[j];
1529  for ( int i = 0; i < 5; i++ ) in >> p[i][j];
1530  for ( int i = 0; i < 15; i++ ) in >> c[i][j];
1531  in >> chi2[j];
1532  in >> ndf[j];
1533  }
1534  t.fX.load( x );
1535  t.fSignCosPhi.load( s );
1536  for ( int i = 0; i < 5; i++ ) t.fP[i].load( p[i] );
1537  for ( int i = 0; i < 5; i++ ) t.fC[i].load( c[i] );
1538  t.fChi2.load( chi2 );
1539  t.fNDF.load( ndf );
1540  return in;
1541 }
Int_t i
Definition: run_full.C:25
TLorentzVector s
Definition: Pnd2DStar.C:50
Double_t p
Definition: anasim.C:58
Double_t x
friend class PndFTSCATrackParam
friend

Definition at line 1364 of file PndFTSCATrackParamVector.h.

Member Data Documentation

float_v PndFTSCATrackParamVector::fAlpha
private
float_v PndFTSCATrackParamVector::fC[15]
private
float_v PndFTSCATrackParamVector::fChi2
private
int_v PndFTSCATrackParamVector::fNDF
private
float_v PndFTSCATrackParamVector::fP[5]
private
float_v PndFTSCATrackParamVector::fSignCosPhi
private
float_v PndFTSCATrackParamVector::fX
private

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