PandaRoot
PndFTSCATrackParam Class Reference

#include <PndFTSCATrackParam.h>

Public Member Functions

 PndFTSCATrackParam ()
 
 PndFTSCATrackParam (const TrackParamVector &v, int i)
 
float X0 () const
 
float X1 () const
 
float X2 () const
 
float X () const
 
float Y () const
 
float Z () const
 
float SinPhi () const
 
float DzDs () const
 
float QPt () const
 
float QP () const
 
float QMomentum () const
 
float SignCosPhi () const
 
float Chi2 () const
 
int NDF () const
 
float Err2X1 () const
 
float Err2X2 () const
 
float Err2Y () const
 
float Err2Z () const
 
float Err2SinPhi () const
 
float Err2DzDs () const
 
float Err2QPt () const
 
float Angle () const
 
float Kappa (float Bz) const
 
float CosPhi () const
 
float Err2QMomentum () const
 
const float * Par () const
 
const float * Cov () const
 
bool GetXYZPxPyPzQ (float &x, float &y, float &z, float &px, float &py, float &pz, int &q, float cov[21]) const
 
void SetSinPhi (float v)
 
void GetDCAPoint (float x, float y, float z, float &px, float &py, float &pz, float Bz) const
 
bool TransportToX0 (float x, float Bz, float maxSinPhi=.999)
 
bool Rotate (float alpha, float maxSinPhi=.999)
 
bool Transport (const FTSCAHit &hit, const PndFTSCAParam &param)
 
bool IsValid () const
 
void SetAsInvalid ()
 

Friends

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

Detailed Description

Definition at line 224 of file PndFTSCATrackParam.h.

Constructor & Destructor Documentation

◆ PndFTSCATrackParam() [1/2]

PndFTSCATrackParam::PndFTSCATrackParam ( )
inline

Definition at line 229 of file PndFTSCATrackParam.h.

229 { Reset(); }

◆ PndFTSCATrackParam() [2/2]

PndFTSCATrackParam::PndFTSCATrackParam ( const TrackParamVector v,
int  i 
)
inline

Definition at line 230 of file PndFTSCATrackParam.h.

References PndCATrackParamVector::Cov(), i, and PndCATrackParamVector::Par().

230  : fX(v.X()[i]), fSignCosPhi(v.SignCosPhi()[i]), fChi2(v.Chi2()[i]), fNDF(v.NDF()[i]), fAlpha(v.Angle()[i])
231  {
232  for (int j = 0; j < 5; ++j)
233  fP[j] = v.Par()[j][i];
234  for (int j = 0; j < 15; ++j)
235  fC[j] = v.Cov()[j][i];
236  }
unsigned int i
Definition: P4_F32vec4.h:21
const float_v & Cov(int i) const
const float_v & Par(int i) const

Member Function Documentation

◆ Angle()

float PndFTSCATrackParam::Angle ( ) const
inline

Definition at line 264 of file PndFTSCATrackParam.h.

264 { return fAlpha; }

◆ Chi2()

float PndFTSCATrackParam::Chi2 ( ) const
inline

Definition at line 252 of file PndFTSCATrackParam.h.

Referenced by FTSCATrack::Fit(), and FTSCATrack::Fit2Times().

252 { return fChi2; }

◆ CosPhi()

float PndFTSCATrackParam::CosPhi ( ) const
inline

Definition at line 267 of file PndFTSCATrackParam.h.

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

Referenced by TransportToX0().

267 { return fSignCosPhi * CAMath::Sqrt(1 - SinPhi() * SinPhi()); }
static T Sqrt(const T &x)
Definition: PndCAMath.h:45

◆ Cov()

const float* PndFTSCATrackParam::Cov ( ) const
inline

Definition at line 272 of file PndFTSCATrackParam.h.

References GetXYZPxPyPzQ(), and pz.

272 { return fC; }

◆ DzDs()

float PndFTSCATrackParam::DzDs ( ) const
inline

Definition at line 246 of file PndFTSCATrackParam.h.

Referenced by QP(), and TransportToX0().

246 { return fP[3]; }

◆ Err2DzDs()

float PndFTSCATrackParam::Err2DzDs ( ) const
inline

Definition at line 261 of file PndFTSCATrackParam.h.

261 { return fC[9]; }

◆ Err2QMomentum()

float PndFTSCATrackParam::Err2QMomentum ( ) const
inline

Definition at line 269 of file PndFTSCATrackParam.h.

Referenced by FTSCANPlet::QMomentumErr(), and FTSCANPlet::QMomentumErr2().

269 { return fC[14]; }

◆ Err2QPt()

float PndFTSCATrackParam::Err2QPt ( ) const
inline

Definition at line 262 of file PndFTSCATrackParam.h.

262 { return fC[14]; }

◆ Err2SinPhi()

float PndFTSCATrackParam::Err2SinPhi ( ) const
inline

Definition at line 260 of file PndFTSCATrackParam.h.

260 { return fC[5]; }

◆ Err2X1()

float PndFTSCATrackParam::Err2X1 ( ) const
inline

Definition at line 255 of file PndFTSCATrackParam.h.

References Err2Y().

255 { return Err2Y(); }

◆ Err2X2()

float PndFTSCATrackParam::Err2X2 ( ) const
inline

Definition at line 256 of file PndFTSCATrackParam.h.

References Err2Z().

256 { return Err2Z(); }

◆ Err2Y()

float PndFTSCATrackParam::Err2Y ( ) const
inline

Definition at line 258 of file PndFTSCATrackParam.h.

Referenced by Err2X1().

258 { return fC[0]; }

◆ Err2Z()

float PndFTSCATrackParam::Err2Z ( ) const
inline

Definition at line 259 of file PndFTSCATrackParam.h.

Referenced by Err2X2().

259 { return fC[2]; }

◆ GetDCAPoint()

void PndFTSCATrackParam::GetDCAPoint ( float  x,
float  y,
float  z,
float &  px,
float &  py,
float &  pz,
float  Bz 
) const

Referenced by SetSinPhi().

◆ GetXYZPxPyPzQ()

bool PndFTSCATrackParam::GetXYZPxPyPzQ ( float &  x,
float &  y,
float &  z,
float &  px,
float &  py,
float &  pz,
int &  q,
float  cov[21] 
) const

Referenced by Cov().

◆ IsValid()

bool PndFTSCATrackParam::IsValid ( ) const
inline

Definition at line 286 of file PndFTSCATrackParam.h.

286 { return fChi2 != -1; }

◆ Kappa()

float PndFTSCATrackParam::Kappa ( float  Bz) const
inline

Definition at line 266 of file PndFTSCATrackParam.h.

266 { return fP[4] * Bz; }

◆ NDF()

int PndFTSCATrackParam::NDF ( ) const
inline

Definition at line 253 of file PndFTSCATrackParam.h.

Referenced by FTSCATrack::Fit(), and FTSCATrack::Fit2Times().

253 { return fNDF; }

◆ Par()

const float* PndFTSCATrackParam::Par ( ) const
inline

Definition at line 271 of file PndFTSCATrackParam.h.

271 { return fP; }

◆ QMomentum()

float PndFTSCATrackParam::QMomentum ( ) const
inline

Definition at line 249 of file PndFTSCATrackParam.h.

References QPt().

Referenced by FTSCANPlet::QMomentum().

249 { return QPt(); }

◆ QP()

float PndFTSCATrackParam::QP ( ) const
inline

Definition at line 248 of file PndFTSCATrackParam.h.

References DzDs(), QPt(), and sqrt().

Referenced by FTSCATrack::Fit().

248 { return QPt() / sqrt(DzDs() * DzDs() + 1); }
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:28

◆ QPt()

float PndFTSCATrackParam::QPt ( ) const
inline

Definition at line 247 of file PndFTSCATrackParam.h.

Referenced by QMomentum(), QP(), and TransportToX0().

247 { return fP[4]; }

◆ Rotate()

bool PndFTSCATrackParam::Rotate ( float  alpha,
float  maxSinPhi = .999 
)

Referenced by SetSinPhi().

◆ SetAsInvalid()

void PndFTSCATrackParam::SetAsInvalid ( )
inline

Definition at line 287 of file PndFTSCATrackParam.h.

287 { fChi2 = -1; }

◆ SetSinPhi()

void PndFTSCATrackParam::SetSinPhi ( float  v)
inline

Definition at line 276 of file PndFTSCATrackParam.h.

References alpha, GetDCAPoint(), Rotate(), Transport(), TransportToX0(), and v.

276 { fP[2] = v; }
__m128 v
Definition: P4_F32vec4.h:3

◆ SignCosPhi()

float PndFTSCATrackParam::SignCosPhi ( ) const
inline

Definition at line 251 of file PndFTSCATrackParam.h.

251 { return fSignCosPhi; }

◆ SinPhi()

float PndFTSCATrackParam::SinPhi ( ) const
inline

Definition at line 245 of file PndFTSCATrackParam.h.

Referenced by CosPhi(), and TransportToX0().

245 { return fP[2]; }

◆ Transport()

bool PndFTSCATrackParam::Transport ( const FTSCAHit hit,
const PndFTSCAParam param 
)

Referenced by SetSinPhi().

◆ TransportToX0()

bool PndFTSCATrackParam::TransportToX0 ( float  x,
float  Bz,
float  maxSinPhi = .999 
)
inline

Definition at line 313 of file PndFTSCATrackParam.h.

References CAMath::Abs(), CAMath::ASin(), CosPhi(), DzDs(), f, QPt(), SinPhi(), CAMath::Sqrt(), X(), Y(), and Z().

Referenced by SetSinPhi().

314 {
315  //* Transport the track parameters to X=x, using linearization at t0, and the field value Bz
316  //* maxSinPhi is the max. allowed value for |t0.SinPhi()|
317  //* linearisation of trajectory t0 is also transported to X=x,
318  //* returns 1 if OK
319  //*
320 
321  const float ex = CosPhi();
322  const float ey = SinPhi();
323  const float k = QPt() * Bz;
324  const float dx = x - X();
325 
326  const float ey1 = k * dx + ey;
327 
328  // check for intersection with X=x
329 
330  if (CAMath::Abs(ey1) > maxSinPhi)
331  return 0;
332 
333  float ex1 = CAMath::Sqrt(1.f - ey1 * ey1);
334  if (ex < 0)
335  ex1 = -ex1;
336 
337  const float dx2 = dx * dx;
338  const float ss = ey + ey1;
339  const float cc = ex + ex1;
340 
341  if ((CAMath::Abs(cc) < 1.e-4 || CAMath::Abs(ex) < 1.e-4 || CAMath::Abs(ex1) < 1.e-4))
342  return 0;
343 
344  const float cci = 1.f / cc;
345  const float exi = 1.f / ex;
346  const float ex1i = 1.f / ex1;
347 
348  const float tg = ss * cci; // tan((phi1+phi)/2)
349 
350  const float dy = dx * tg;
351  float dl = dx * CAMath::Sqrt(1.f + tg * tg);
352 
353  if (cc < 0)
354  dl = -dl;
355  float dSin = dl * k * 0.5;
356  if (dSin > 1.f)
357  dSin = 1.f;
358  if (dSin < -1.f)
359  dSin = -1.f;
360  const float dS = (CAMath::Abs(k) > 1.e-4) ? (2 * CAMath::ASin(dSin) / k) : dl;
361  const float dz = dS * DzDs();
362 
363  // float H0[5] = { 1,0, h2, 0, h4 };
364  // float H1[5] = { 0, 1, 0, dS, 0 };
365  // float H2[5] = { 0, 0, 1, 0, dxBz };
366  // float H3[5] = { 0, 0, 0, 1, 0 };
367  // float H4[5] = { 0, 0, 0, 0, 1 };
368 
369  const float h2 = dx * (1.f + ey * ey1 + ex * ex1) * exi * ex1i * cci;
370  const float h4 = dx2 * (cc + ss * ey1 * ex1i) * cci * cci * Bz;
371  const float dxBz = dx * Bz;
372 
373  fX = X() + dx;
374  fP[0] = Y() + dy;
375  fP[1] = Z() + dz;
376  fP[2] = ey1;
377 
378 #if 1
379  const float c00 = fC[0];
380  const float c10 = fC[1];
381  const float c11 = fC[2];
382  const float c20 = fC[3];
383  const float c21 = fC[4];
384  const float c22 = fC[5];
385  const float c30 = fC[6];
386  const float c31 = fC[7];
387  const float c32 = fC[8];
388  const float c33 = fC[9];
389  const float c40 = fC[10];
390  const float c41 = fC[11];
391  const float c42 = fC[12];
392  const float c43 = fC[13];
393  const float c44 = fC[14];
394 
395  fC[0] = (c00 + h2 * h2 * c22 + h4 * h4 * c44 + 2.f * (h2 * c20 + h4 * c40 + h2 * h4 * c42));
396 
397  fC[1] = c10 + h2 * c21 + h4 * c41 + dS * (c30 + h2 * c32 + h4 * c43);
398  fC[2] = c11 + 2.f * dS * c31 + dS * dS * c33;
399 
400  fC[3] = c20 + h2 * c22 + h4 * c42 + dxBz * (c40 + h2 * c42 + h4 * c44);
401  fC[4] = c21 + dS * c32 + dxBz * (c41 + dS * c43);
402  fC[5] = c22 + 2.f * dxBz * c42 + dxBz * dxBz * c44;
403 
404  fC[6] = c30 + h2 * c32 + h4 * c43;
405  fC[7] = c31 + dS * c33;
406  fC[8] = c32 + dxBz * c43;
407  fC[9] = c33;
408 
409  fC[10] = c40 + h2 * c42 + h4 * c44;
410  fC[11] = c41 + dS * c43;
411  fC[12] = c42 + dxBz * c44;
412  fC[13] = c43;
413  fC[14] = c44;
414 #else
415  fC[0] = (fC[0] + h2 * h2 * fC[5] + h4 * h4 * fC[14] + 2 * (h2 * fC[3] + h4 * fC[10] + h2 * h4 * fC[12]));
416 
417  fC[1] = fC[1] + h2 * fC[4] + h4 * fC[11] + dS * (fC[6] + h2 * fC[8] + h4 * fC[13]);
418  fC[2] = fC[2] + 2 * dS * fC[7] + dS * dS * fC[9];
419 
420  fC[3] = fC[3] + h2 * fC[5] + h4 * fC[12] + dxBz * (fC[10] + h2 * fC[12] + h4 * fC[14]);
421  fC[4] = fC[4] + dS * fC[8] + dxBz * (fC[11] + dS * fC[13]);
422  fC[5] = fC[5] + 2 * dxBz * fC[12] + dxBz * dxBz * fC[14];
423 
424  fC[6] = fC[6] + h2 * fC[8] + h4 * fC[13];
425  fC[7] = fC[7] + dS * fC[9];
426  fC[8] = fC[8] + dxBz * fC[13];
427  fC[9] = fC[9];
428 
429  fC[10] = fC[10] + h2 * fC[12] + h4 * fC[14];
430  fC[11] = fC[11] + dS * fC[13];
431  fC[12] = fC[12] + dxBz * fC[14];
432  fC[13] = fC[13];
433  fC[14] = fC[14];
434 #endif
435 
436  return 1;
437 }
static T ASin(const T &x)
static T Sqrt(const T &x)
Definition: PndCAMath.h:45
static T Abs(const T &x)
Definition: PndCAMath.h:56
float f
Definition: P4_F32vec4.h:20

◆ X()

float PndFTSCATrackParam::X ( ) const
inline

Definition at line 242 of file PndFTSCATrackParam.h.

Referenced by TransportToX0().

242 { return fX; }

◆ X0()

float PndFTSCATrackParam::X0 ( ) const
inline

Definition at line 238 of file PndFTSCATrackParam.h.

238 { return fX; }

◆ X1()

float PndFTSCATrackParam::X1 ( ) const
inline

Definition at line 239 of file PndFTSCATrackParam.h.

References Y().

239 { return Y(); }

◆ X2()

float PndFTSCATrackParam::X2 ( ) const
inline

Definition at line 240 of file PndFTSCATrackParam.h.

References Z().

240 { return Z(); }

◆ Y()

float PndFTSCATrackParam::Y ( ) const
inline

Definition at line 243 of file PndFTSCATrackParam.h.

Referenced by TransportToX0(), and X1().

243 { return fP[0]; }

◆ Z()

float PndFTSCATrackParam::Z ( ) const
inline

Definition at line 244 of file PndFTSCATrackParam.h.

Referenced by TransportToX0(), and X2().

244 { return fP[1]; }

Friends And Related Function Documentation

◆ operator<<

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

◆ operator>>

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

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