WTrackParameter Class Reference

#include <WTrackParameter.h>

List of all members.

Public Member Functions

 WTrackParameter ()
 ~WTrackParameter ()
 WTrackParameter (const WTrackParameter &wtrk)
WTrackParameteroperator= (const WTrackParameter &wtrk)
 WTrackParameter (const int charge, const HepLorentzVector &p, const HepPoint3D &x, const double err[])
 WTrackParameter (const int charge, const HepLorentzVector &p, const HepPoint3D &x, const HepSymMatrix &err)
 WTrackParameter (const double mass, const HepVector &helix, const double err[])
 WTrackParameter (const double mass, const HepVector &helix, const HepSymMatrix &err)
 WTrackParameter (const HepPoint3D &x, const HepLorentzVector &p, const double dphi, const double dtheta, const double dE)
 WTrackParameter (const HepLorentzVector &p, const double dphi, const double dtheta, const double dE)
void setW (const HepVector &w)
void setW (const int n, const double w)
void setEw (const HepSymMatrix &Ew)
void setCharge (const int charge)
void setMass (const double mass)
void setType (const int type)
void setVplm (const HepSymMatrix &Vplm)
void setPlmp (const HepVector &plmp)
int type () const
int charge () const
double mass () const
bool IsInvariableMass () const
double phi () const
double Lambda () const
HepVector w () const
HepSymMatrix Ew () const
HepLorentzVector p () const
HepSymMatrix Ep () const
HepPoint3D x () const
HepVector X () const
HepSymMatrix Ex () const
HepVector philambdamass ()
HepSymMatrix Vplm () const
HepVector plmp () const

Private Member Functions

HepMatrix GetCvtMatrix (const double mass, const HepVector &helix)
HepVector CvtH2W (const double mass, const HepVector &helix)

Private Attributes

int m_charge
HepVector m_w
HepSymMatrix m_Ew
int m_type
HepSymMatrix m_Vplm
HepVector m_plmp
double m_mass
bool m_massInvariable


Detailed Description

Definition at line 31 of file WTrackParameter.h.


Constructor & Destructor Documentation

WTrackParameter::WTrackParameter (  ) 

Definition at line 3 of file WTrackParameter.cxx.

References m_charge, m_Ew, m_mass, m_massInvariable, m_type, and m_w.

Referenced by WTrackParameter().

00004 {
00005         m_w      = HepVector(7, 0);
00006         m_Ew     = HepSymMatrix(7, 0);
00007         m_charge = 0;
00008         m_type   = 1;
00009         m_mass   = 0;
00010         m_massInvariable = true;
00011 }

WTrackParameter::~WTrackParameter (  )  [inline]

Definition at line 36 of file WTrackParameter.h.

00036 {;}

WTrackParameter::WTrackParameter ( const WTrackParameter wtrk  ) 

Definition at line 13 of file WTrackParameter.cxx.

References m_charge, m_Ew, m_mass, m_massInvariable, m_plmp, m_type, m_Vplm, and m_w.

00014 {
00015         m_charge = wtrk.m_charge;
00016         m_w      = wtrk.m_w;
00017         m_Ew     = wtrk.m_Ew;
00018         m_type   = wtrk.m_type;
00019         m_plmp   = wtrk.m_plmp;
00020         m_Vplm   = wtrk.m_Vplm;
00021         m_mass   = wtrk.m_mass;
00022         m_massInvariable = wtrk.m_massInvariable;
00023 }

WTrackParameter::WTrackParameter ( const int  charge,
const HepLorentzVector &  p,
const HepPoint3D x,
const double  err[] 
)

Definition at line 53 of file WTrackParameter.cxx.

References ers::error, genRecEmupikp::i, ganga-rec::j, and WTrackParameter().

00054 {
00055         HepSymMatrix error(5, 0);
00056         int k = 0;
00057         for (int i = 0; i < 5; i++) 
00058         {
00059                 for(int j = i; j < 5; j++) 
00060                 {
00061                         error[i][j] = err[k];
00062                         error[j][i] = err[k];
00063                         k++;
00064                 }
00065         }
00066         *this = WTrackParameter(charge, p, x, error);
00067 }

WTrackParameter::WTrackParameter ( const int  charge,
const HepLorentzVector &  p,
const HepPoint3D x,
const HepSymMatrix &  err 
)

Definition at line 69 of file WTrackParameter.cxx.

References cos(), GetCvtMatrix(), genRecEmupikp::i, m_charge, m_Ew, m_mass, m_massInvariable, m_type, m_w, phi0, and sin().

00070 {
00071         m_w      = HepVector(7, 0);
00072         m_Ew     = HepSymMatrix(7, 0);
00073         m_charge = charge;
00074         m_type   = 1;
00075         m_mass   = p.m();
00076         m_massInvariable = true;
00077         for (int i = 0; i < 4; i++)
00078                 m_w[i] = p[i];
00079         for (int i = 0; i < 3; i++)
00080                 m_w[i+4] = x[i];
00081 
00082         HepVector helix(5, 0);
00083         double phi0 = atan2(-p[0], p[1]);
00084         if (cos(phi0) != 0)
00085                 helix[0] = x[0] / cos(phi0);
00086         else
00087                 helix[0] = x[1] / sin(phi0);
00088         helix[1] = phi0;
00089         helix[2] = charge/p.perp();
00090         helix[3] = x[3];
00091         helix[4] = p[3]/p.perp();
00092         HepMatrix dWdA(7, 5, 0);
00093         dWdA = GetCvtMatrix(p.m(), helix);
00094         m_Ew = err.similarity(dWdA);
00095 }

WTrackParameter::WTrackParameter ( const double  mass,
const HepVector &  helix,
const double  err[] 
)

Definition at line 97 of file WTrackParameter.cxx.

References ers::error, genRecEmupikp::i, ganga-rec::j, and WTrackParameter().

00098 { 
00099         HepSymMatrix error(5, 0);
00100         int k = 0;
00101         for(int i = 0; i < 5; i++) 
00102         {
00103                 for(int j = i; j < 5; j++) 
00104                 {
00105                         error[i][j] = err[k];
00106                         error[j][i] = err[k];
00107                         k++;
00108                 }
00109         }
00110         *this = WTrackParameter(mass, hel, error);
00111 }

WTrackParameter::WTrackParameter ( const double  mass,
const HepVector &  helix,
const HepSymMatrix &  err 
)

Definition at line 114 of file WTrackParameter.cxx.

References CvtH2W(), GetCvtMatrix(), m_charge, m_Ew, m_mass, m_massInvariable, m_type, and m_w.

00115 {
00116         m_type   = 1;
00117         m_mass   = mass;
00118         m_w      = HepVector(7, 0);
00119         m_Ew     = HepSymMatrix(7, 0);
00120         HepMatrix dWdA(7, 5, 0);
00121         dWdA     = GetCvtMatrix(mass, hel);
00122         m_Ew     = err.similarity(dWdA);
00123         m_w      = CvtH2W(mass, hel);
00124         m_charge = hel[2] > 0 ? +1 : -1;
00125         m_massInvariable = true;
00126 }

WTrackParameter::WTrackParameter ( const HepPoint3D x,
const HepLorentzVector &  p,
const double  dphi,
const double  dtheta,
const double  dE 
)

Definition at line 135 of file WTrackParameter.cxx.

References genRecEmupikp::i, m_charge, m_Ew, m_mass, m_massInvariable, m_type, and m_w.

00136 {
00137         m_w      = HepVector(7, 0);
00138         m_Ew     = HepSymMatrix(7, 0);
00139         m_type   = 2;
00140         m_mass   = p.m();
00141         m_charge = 0;
00142         m_massInvariable = true;
00143         for (int i = 0; i< 4; i++)
00144                 m_w[i] = p[i];
00145         for (int i = 0; i< 3; i++)
00146                 m_w[i+4] = x[i];
00147         
00148         HepMatrix dwda(7, 3, 0);
00149         dwda[0][0] = -p.py();
00150         dwda[1][0] = p.px();
00151         dwda[4][0] = -x.y();
00152         dwda[5][0] = x.x();
00153 
00154         dwda[0][1] = p.px()*p.pz()/p.perp();
00155         dwda[1][1] = p.py()*p.pz()/p.perp();
00156         dwda[2][1] = -p.perp();
00157         dwda[6][1] = -x.distance2()/sqrt(x.x()*x.x()+x.y()*x.y());
00158 
00159         dwda[0][2] = p.px()/p.rho();
00160         dwda[1][2] = p.py()/p.rho();
00161         dwda[2][2] = p.pz()/p.rho();
00162         dwda[3][2] = p.rho()/p.e();
00163         // dwda[4][2] = x.x()/(x.x()*x.x()+x.y()*x.y());
00164         // dwda[5][2] = x.y()/(x.x()*x.x()+x.y()*x.y());
00165         // dwda[6][2] = x.z()/(x.x()*x.x()+x.y()*x.y());
00166 
00167         HepSymMatrix emcmea(3, 0);
00168         emcmea[0][0] = dphi * dphi;
00169         emcmea[1][1] = dthe * dthe;
00170         emcmea[2][2] = dE * dE;
00171         m_Ew = emcmea.similarity(dwda);
00172 }

WTrackParameter::WTrackParameter ( const HepLorentzVector &  p,
const double  dphi,
const double  dtheta,
const double  dE 
)

Definition at line 128 of file WTrackParameter.cxx.

References WTrackParameter(), and x().

00129 {
00130         HepPoint3D x(0, 0, 0);
00131         *this = WTrackParameter(x, p, dphi, dthe, dE);
00132 }


Member Function Documentation

int WTrackParameter::charge (  )  const [inline]

Definition at line 62 of file WTrackParameter.h.

References m_charge.

Referenced by HTrackParameter::HTrackParameter(), VertexFit::pull(), and VertexConstraints::UpdateConstraints().

00062 {return m_charge;}

HepVector WTrackParameter::CvtH2W ( const double  mass,
const HepVector &  helix 
) [private]

Definition at line 205 of file WTrackParameter.cxx.

References cos(), lambda, phi0, q, sin(), and w().

Referenced by WTrackParameter().

00206 {
00207         double drho    = helix[0];
00208         double phi0    = helix[1];
00209         double kappa   = helix[2];
00210         double dz      = helix[3];
00211         double lambda  = helix[4];
00212         double sinphi0 = sin(phi0);
00213         double cosphi0 = cos(phi0);
00214         int    q       = (kappa>0) ? 1 : (-1);
00215         HepVector w(7, 0);
00216         w[0] = -sinphi0 * q / kappa;
00217         w[1] = cosphi0 * q / kappa;
00218         w[2] = lambda * q / kappa;
00219         w[3] = sqrt( (1+lambda*lambda) / (kappa*kappa) + mass * mass );
00220         w[4] = drho * cosphi0;
00221         w[5] = drho * sinphi0;
00222         w[6] = dz;
00223         return w;
00224 }

HepSymMatrix WTrackParameter::Ep (  )  const [inline]

Definition at line 70 of file WTrackParameter.h.

References m_Ew.

00070 {return m_Ew.sub(1, 4); }

HepSymMatrix WTrackParameter::Ew (  )  const [inline]

Definition at line 68 of file WTrackParameter.h.

References m_Ew.

Referenced by TrackPool::AddMissTrack(), TrackPool::AddTrack(), LambdaReconstruction::execute(), KShortReconstruction::execute(), and HTrackParameter::HTrackParameter().

00068 {return m_Ew;}

HepSymMatrix WTrackParameter::Ex (  )  const [inline]

Definition at line 73 of file WTrackParameter.h.

References m_Ew.

00073 {return m_Ew.sub(5, 7); }

HepMatrix WTrackParameter::GetCvtMatrix ( const double  mass,
const HepVector &  helix 
) [private]

Definition at line 174 of file WTrackParameter.cxx.

References cos(), lambda, phi0, q, and sin().

Referenced by WTrackParameter().

00175 {
00176         HepMatrix m(7, 5, 0);
00177         double drho    = helix[0];
00178         double phi0    = helix[1];
00179         double kappa   = helix[2];
00180         double kappa2  = kappa * kappa;
00181         double kappa3  = kappa * kappa2;
00182         double dz      = helix[3];
00183         double lambda  = helix[4];
00184         double lambda2 = lambda * lambda;
00185         double e       = sqrt( (1+lambda2) / kappa2 + mass * mass );
00186         double sinphi0 = sin(phi0);
00187         double cosphi0 = cos(phi0);
00188         int    q       = (kappa>0) ? 1 : (-1);
00189         m[0][1] = -cosphi0 * q / kappa;
00190         m[0][2] = sinphi0 * q / kappa2;
00191         m[1][1] = -sinphi0 * q / kappa;
00192         m[1][2] = -cosphi0 * q / kappa2;
00193         m[2][2] = -lambda * q / kappa2;
00194         m[2][4] = q / kappa;
00195         m[3][2] = -(1+lambda2) / (kappa3 * e);
00196         m[3][4] = lambda / (kappa2 * e);
00197         m[4][0] = cosphi0;
00198         m[4][1] = -drho * sinphi0;
00199         m[5][0] = sinphi0;
00200         m[5][1] = drho * cosphi0;
00201         m[6][3] = 1;
00202         return m;
00203 }

bool WTrackParameter::IsInvariableMass (  )  const [inline]

Definition at line 64 of file WTrackParameter.h.

References m_massInvariable.

00064 {return m_massInvariable;}

double WTrackParameter::Lambda (  )  const [inline]

Definition at line 66 of file WTrackParameter.h.

References w().

Referenced by TrackPool::AddMissTrack().

00066 {return w()[2]/sqrt(w()[0]*w()[0] + w()[1]*w()[1]);}

double WTrackParameter::mass (  )  const [inline]

Definition at line 63 of file WTrackParameter.h.

References m_mass.

Referenced by VertexFit::pull(), and VertexFit::swimVertex().

00063 {return m_mass;} 

WTrackParameter & WTrackParameter::operator= ( const WTrackParameter wtrk  ) 

Definition at line 25 of file WTrackParameter.cxx.

References m_charge, m_Ew, m_mass, m_massInvariable, m_plmp, m_type, m_Vplm, and m_w.

00026 {
00027         if (this == &wtrk)
00028                 return *this;
00029         m_charge = wtrk.m_charge;
00030         m_type   = wtrk.m_type;
00031         m_w      = wtrk.m_w;
00032         m_Ew     = wtrk.m_Ew;
00033         m_plmp   = wtrk.m_plmp;
00034         m_Vplm   = wtrk.m_Vplm;
00035         m_mass   = wtrk.m_mass;
00036         m_massInvariable = wtrk.m_massInvariable;
00037         return *this;
00038 }

HepLorentzVector WTrackParameter::p (  )  const [inline]

Definition at line 69 of file WTrackParameter.h.

References m_w.

Referenced by LambdaReconstruction::execute(), KShortReconstruction::execute(), DQASelHadron::execute(), DQASelDimu::execute(), DQASelBhabha::execute(), DQAKsKpi::execute(), HTrackParameter::HTrackParameter(), philambdamass(), and VertexConstraints::UpdateConstraints().

00069 {return HepLorentzVector(m_w[0], m_w[1], m_w[2], m_w[3]);}

double WTrackParameter::phi (  )  const [inline]

Definition at line 65 of file WTrackParameter.h.

References w().

00065 {return atan(w()[1]/(w()[0]));}

HepVector WTrackParameter::philambdamass (  ) 

Definition at line 40 of file WTrackParameter.cxx.

References p(), and w().

00041 {
00042         HepVector tmp(4,0);
00043         if (w()[1] >= 0)
00044                 tmp[0] = atan(w()[1]/(w()[0]));
00045         else 
00046                 tmp[0] = atan(w()[1]/(w()[0])) + 3.1415926;
00047         tmp[1] = w()[2]/sqrt(w()[0]*w()[0] + w()[1]*w()[1]);
00048         tmp[2] = p().m();
00049         tmp[3] = sqrt(w()[3]*w()[3] - p().m()*p().m());
00050         return tmp;
00051 }

HepVector WTrackParameter::plmp (  )  const [inline]

Definition at line 76 of file WTrackParameter.h.

References m_plmp.

00076 {return m_plmp;}

void WTrackParameter::setCharge ( const int  charge  )  [inline]

Definition at line 54 of file WTrackParameter.h.

References m_charge.

Referenced by VertexFit::BuildVirtualParticle(), KinematicFit::BuildVirtualParticle(), KalmanKinematicFit::BuildVirtualParticle(), SecondVertexFit::Fit(), VertexFit::pull(), KalmanVertexFit::wTrack(), and HTrackParameter::wTrack().

00054 {m_charge = charge;}

void WTrackParameter::setEw ( const HepSymMatrix &  Ew  )  [inline]

Definition at line 53 of file WTrackParameter.h.

References m_Ew.

Referenced by TrackPool::AddMissTrack(), TrackPool::AddTrack(), VertexFit::BuildVirtualParticle(), KinematicFit::BuildVirtualParticle(), KalmanKinematicFit::BuildVirtualParticle(), KinematicFit::covMatrix(), SecondVertexFit::Fit(), VertexFit::pull(), KinematicFit::pull(), KalmanKinematicFit::pull(), VertexFit::swimVertex(), KinematicFit::upCovmtx(), KalmanVertexFit::wTrack(), and HTrackParameter::wTrack().

00053 {m_Ew = Ew;}

void WTrackParameter::setMass ( const double  mass  )  [inline]

Definition at line 55 of file WTrackParameter.h.

References m_mass.

Referenced by TrackPool::AddMissTrack(), KinematicFit::BuildVirtualParticle(), and KalmanKinematicFit::BuildVirtualParticle().

00055 {m_mass = mass;}

void WTrackParameter::setPlmp ( const HepVector &  plmp  )  [inline]

Definition at line 58 of file WTrackParameter.h.

References m_plmp.

Referenced by TrackPool::AddMissTrack(), and TrackPool::AddTrack().

00058 {m_plmp = plmp;}

void WTrackParameter::setType ( const int  type  )  [inline]

Definition at line 56 of file WTrackParameter.h.

References m_type.

Referenced by TrackPool::AddMissTrack().

00056 {m_type = type;}

void WTrackParameter::setVplm ( const HepSymMatrix &  Vplm  )  [inline]

Definition at line 57 of file WTrackParameter.h.

References m_Vplm.

Referenced by TrackPool::AddMissTrack(), and TrackPool::AddTrack().

00057 {m_Vplm = Vplm;}

void WTrackParameter::setW ( const int  n,
const double  w 
) [inline]

Definition at line 52 of file WTrackParameter.h.

References m_w.

00052 {m_w[n] = w;}

void WTrackParameter::setW ( const HepVector &  w  )  [inline]

Definition at line 51 of file WTrackParameter.h.

References m_mass, and m_w.

Referenced by TrackPool::AddMissTrack(), TrackPool::AddTrack(), VertexFit::BuildVirtualParticle(), KinematicFit::BuildVirtualParticle(), KalmanKinematicFit::BuildVirtualParticle(), SecondVertexFit::Fit(), VertexFit::pull(), KinematicFit::pull(), KalmanKinematicFit::pull(), VertexFit::swimVertex(), KinematicFit::upCovmtx(), KalmanVertexFit::wTrack(), and HTrackParameter::wTrack().

00051 {m_w = w; m_mass = sqrt(w[3]*w[3] - w[2]*w[2] - w[1]*w[1] - w[0]*w[0]);}

int WTrackParameter::type (  )  const [inline]

Definition at line 61 of file WTrackParameter.h.

References m_type.

00061 {return m_type;}

HepSymMatrix WTrackParameter::Vplm (  )  const [inline]

Definition at line 75 of file WTrackParameter.h.

References m_Vplm.

00075 {return m_Vplm;}

HepVector WTrackParameter::w (  )  const [inline]

Definition at line 67 of file WTrackParameter.h.

References m_w.

Referenced by CvtH2W(), LambdaReconstruction::execute(), KShortReconstruction::execute(), HTrackParameter::HTrackParameter(), Lambda(), Pipipi0::MTotal(), Pipi::MTotal(), Kpipi0pi0::MTotal(), Kpipi0::MTotal(), Kpi::MTotal(), Kkpipi::MTotal(), Kkpi0::MTotal(), Kk::MTotal(), K3pipi0::MTotal(), K3pi::MTotal(), K0pipipi0::MTotal(), K0pipi::MTotal(), K0pi0::MTotal(), K0kpi::MTotal(), K0kk::MTotal(), phi(), and philambdamass().

00067 {return m_w;}

HepVector WTrackParameter::X (  )  const [inline]

Definition at line 72 of file WTrackParameter.h.

References m_w.

Referenced by VertexConstraints::UpdateConstraints().

00072 {return m_w.sub(5, 7); }

HepPoint3D WTrackParameter::x (  )  const [inline]

Definition at line 71 of file WTrackParameter.h.

References m_w.

Referenced by VertexConstraints::UpdateConstraints(), and WTrackParameter().

00071 {return HepPoint3D(m_w[4], m_w[5], m_w[6]);}


Member Data Documentation

int WTrackParameter::m_charge [private]

Definition at line 80 of file WTrackParameter.h.

Referenced by charge(), operator=(), setCharge(), and WTrackParameter().

HepSymMatrix WTrackParameter::m_Ew [private]

Definition at line 82 of file WTrackParameter.h.

Referenced by Ep(), Ew(), Ex(), operator=(), setEw(), and WTrackParameter().

double WTrackParameter::m_mass [private]

Definition at line 86 of file WTrackParameter.h.

Referenced by mass(), operator=(), setMass(), setW(), and WTrackParameter().

bool WTrackParameter::m_massInvariable [private]

Definition at line 87 of file WTrackParameter.h.

Referenced by IsInvariableMass(), operator=(), and WTrackParameter().

HepVector WTrackParameter::m_plmp [private]

Definition at line 85 of file WTrackParameter.h.

Referenced by operator=(), plmp(), setPlmp(), and WTrackParameter().

int WTrackParameter::m_type [private]

Definition at line 83 of file WTrackParameter.h.

Referenced by operator=(), setType(), type(), and WTrackParameter().

HepSymMatrix WTrackParameter::m_Vplm [private]

Definition at line 84 of file WTrackParameter.h.

Referenced by operator=(), setVplm(), Vplm(), and WTrackParameter().

HepVector WTrackParameter::m_w [private]

Definition at line 81 of file WTrackParameter.h.

Referenced by operator=(), p(), setW(), w(), WTrackParameter(), X(), and x().


Generated on Tue Nov 29 23:36:21 2016 for BOSS_7.0.2 by  doxygen 1.4.7