Helix Class Reference

Helix parameter class. More...

#include <Helix.h>

List of all members.

Public Member Functions

 Helix (const HepPoint3D &pivot, const HepVector &a, const HepSymMatrix &Ea)
 Constructor with pivot, helix parameter a, and its error matrix.
 Helix (const HepPoint3D &pivot, const HepVector &a)
 Constructor without error matrix.
 Helix (const HepPoint3D &position, const Hep3Vector &momentum, double charge)
 Constructor with position, momentum, and charge.
virtual ~Helix ()
 Destructor.
const HepPoint3Dcenter (void) const
 returns position of helix center(z = 0.);
const HepPoint3Dpivot (void) const
 returns pivot position.
double radius (void) const
 returns radious of helix.
HepPoint3D x (double dPhi=0.) const
 returns position after rotating angle dPhi in phi direction.
double * x (double dPhi, double p[3]) const
HepPoint3D x (double dPhi, HepSymMatrix &Ex) const
 returns position and convariance matrix(Ex) after rotation.
Hep3Vector direction (double dPhi=0.) const
 returns direction vector after rotating angle dPhi in phi direction.
Hep3Vector momentum (double dPhi=0.) const
 returns momentum vector after rotating angle dPhi in phi direction.
Hep3Vector momentum (double dPhi, HepSymMatrix &Em) const
 returns momentum vector after rotating angle dPhi in phi direction.
HepLorentzVector momentum (double dPhi, double mass) const
 returns 4momentum vector after rotating angle dPhi in phi direction.
HepLorentzVector momentum (double dPhi, double mass, HepSymMatrix &Em) const
 returns 4momentum vector after rotating angle dPhi in phi direction.
HepLorentzVector momentum (double dPhi, double mass, HepPoint3D &x, HepSymMatrix &Emx) const
 returns 4momentum vector after rotating angle dPhi in phi direction.
double dr (void) const
 returns an element of parameters.
double phi0 (void) const
double kappa (void) const
double dz (void) const
double tanl (void) const
double curv (void) const
double sinPhi0 (void) const
double cosPhi0 (void) const
const HepVector & a (void) const
 returns helix parameters.
const HepSymMatrix & Ea (void) const
 returns error matrix.
double pt (void) const
double cosTheta (void) const
const HepVector & a (const HepVector &newA)
 sets helix parameters.
const HepSymMatrix & Ea (const HepSymMatrix &newdA)
 sets helix paramters and error matrix.
const HepPoint3Dpivot (const HepPoint3D &newPivot)
 sets pivot position.
void set (const HepPoint3D &pivot, const HepVector &a, const HepSymMatrix &Ea)
 sets helix pivot position, parameters, and error matrix.
void ignoreErrorMatrix (void)
 unsets error matrix. Error calculations will be ignored after this function call until an error matrix be set again. 0 matrix will be return as a return value for error matrix when you call functions which returns an error matrix.
double bFieldZ (double)
 sets/returns z componet of the magnetic field.
double bFieldZ (void) const
Helixoperator= (const Helix &)
 Copy operator.
HepMatrix delApDelA (const HepVector &ap) const
HepMatrix delXDelA (double phi) const
HepMatrix delMDelA (double phi) const
HepMatrix del4MDelA (double phi, double mass) const
HepMatrix del4MXDelA (double phi, double mass) const

Static Public Attributes

static const double ConstantAlpha = 333.564095
 Constant alpha for uniform field.

Protected Attributes

IMagneticFieldSvcm_pmgnIMF
double m_bField
double m_alpha

Private Member Functions

void updateCache (void)

Private Attributes

HepPoint3D m_pivot
HepVector m_a
HepSymMatrix m_Ea
bool m_matrixValid
HepPoint3D m_center
double m_cp
double m_sp
double m_pt
double m_r
double m_ac [5]


Detailed Description

Helix parameter class.

Definition at line 53 of file Helix.h.


Constructor & Destructor Documentation

Helix::Helix ( const HepPoint3D pivot,
const HepVector &  a,
const HepSymMatrix &  Ea 
)

Constructor with pivot, helix parameter a, and its error matrix.

Definition at line 47 of file Helix.cxx.

References IMagneticFieldSvc::getReferField(), m_alpha, m_bField, m_pmgnIMF, and updateCache().

00050   : //m_bField(-10.0),
00051   //m_alpha(-333.564095),
00052   m_pivot(pivot),
00053   m_a(a),
00054   m_matrixValid(true),
00055   m_Ea(Ea) {
00056   StatusCode scmgn = Gaudi::svcLocator()->service("MagneticFieldSvc",m_pmgnIMF); 
00057   if(scmgn!=StatusCode::SUCCESS) { 
00058     std::cout<< "Unable to open Magnetic field service"<<std::endl;
00059   }
00060   m_bField = -10000*(m_pmgnIMF->getReferField());
00061   m_alpha = 10000. / 2.99792458 / m_bField;
00062   // m_alpha = 10000. / 2.99792458 / m_bField;
00063   // m_alpha = 333.564095;
00064   updateCache();
00065 }

Helix::Helix ( const HepPoint3D pivot,
const HepVector &  a 
)

Constructor without error matrix.

Definition at line 67 of file Helix.cxx.

References IMagneticFieldSvc::getReferField(), m_alpha, m_bField, m_pmgnIMF, and updateCache().

00069   : //m_bField(-10.0),
00070   //m_alpha(-333.564095),
00071   m_pivot(pivot),
00072   m_a(a),
00073   m_matrixValid(false),
00074   m_Ea(HepSymMatrix(5,0)) {
00075   StatusCode scmgn = Gaudi::svcLocator()->service("MagneticFieldSvc",m_pmgnIMF); 
00076   if(scmgn!=StatusCode::SUCCESS) { 
00077     // log << MSG::ERROR << "Unable to open Magnetic field service"<<endreq; 
00078     std::cout<< "Unable to open Magnetic field service"<<std::endl;
00079   }
00080   m_bField = -10000*(m_pmgnIMF->getReferField());
00081   m_alpha = 10000. / 2.99792458 / m_bField;
00082     // m_alpha = 333.564095;
00083  //cout<<"MdcFastTrakAlg:: bField,alpha: "<<m_bField<<" , "<<m_alpha<<endl;
00084     updateCache();
00085 }

Helix::Helix ( const HepPoint3D position,
const Hep3Vector &  momentum,
double  charge 
)

Constructor with position, momentum, and charge.

Definition at line 87 of file Helix.cxx.

References DBL_MAX, IMagneticFieldSvc::getReferField(), m_a, m_alpha, m_bField, M_PI2, M_PI4, m_pmgnIMF, and updateCache().

00090   : //m_bField(-10.0),
00091   //m_alpha(-333.564095),
00092   m_pivot(position),
00093   m_a(HepVector(5,0)),
00094   m_matrixValid(false),
00095   m_Ea(HepSymMatrix(5,0)) {
00096   StatusCode scmgn = Gaudi::svcLocator()->service("MagneticFieldSvc",m_pmgnIMF); 
00097   if(scmgn!=StatusCode::SUCCESS) { 
00098     // log << MSG::ERROR << "Unable to open Magnetic field service"<<endreq; 
00099     std::cout<< "Unable to open Magnetic field service"<<std::endl;
00100   }
00101   m_bField = -10000*(m_pmgnIMF->getReferField());
00102   m_alpha = 10000. / 2.99792458 / m_bField;
00103   
00104   m_a[0] = 0.;
00105     m_a[1] = fmod(atan2(- momentum.x(), momentum.y())
00106                   + M_PI4, M_PI2);
00107     m_a[3] = 0.;
00108     double perp(momentum.perp());
00109     if (perp != 0.0) {
00110         m_a[2] = charge / perp;
00111         m_a[4] = momentum.z() / perp; 
00112     }
00113     else {
00114         m_a[2] = charge * (DBL_MAX);
00115         if (momentum.z() >= 0) {
00116             m_a[4] = (DBL_MAX);
00117         } else {
00118             m_a[4] = -(DBL_MAX);
00119         }
00120     }
00121     // m_alpha = 333.564095;
00122     updateCache();
00123 }

Helix::~Helix (  )  [virtual]

Destructor.

Definition at line 125 of file Helix.cxx.

00125               {
00126 }


Member Function Documentation

const HepVector & Helix::a ( const HepVector &  newA  )  [inline]

sets helix parameters.

Definition at line 275 of file Helix.h.

References m_a, and updateCache().

00275                             {
00276     m_a = i;
00277     updateCache();
00278     return m_a;
00279 }

const HepVector & Helix::a ( void   )  const [inline]

returns helix parameters.

Definition at line 263 of file Helix.h.

References m_a.

Referenced by TTrack::approach2D(), TBuilderCurl::buildStereo(), TBuilderCosmic::buildStereo(), TBuilder0::buildStereo(), TBuilder::buildStereo(), TBuilder0::buildStereo0(), TBuilderCurl::buildStereoMC(), TCurlFinder::distance(), TTrack::dump(), TTrack::dxda(), TTrack::dxda2D(), Emc_helix::Emc_Get(), EsTimeAlg::execute(), TTrack::fit2D(), MdcHoughFinder::GetMcInfo(), HoughValidUpdate::GetMcInfo(), TTrack::HelCyl(), HelixHasNan(), FTFinder::makeTds(), TRunge::pivot(), TTrack::pt(), TTrack::ptot(), TTrack::pz(), TofFz_helix::TofFz_Get(), TrackKinematics(), and TTrack::TTrack().

00263                    {
00264     return m_a;
00265 }

double Helix::bFieldZ ( void   )  const [inline]

Definition at line 298 of file Helix.h.

References m_bField.

00298                          {
00299     return m_bField;
00300 }

double Helix::bFieldZ ( double   )  [inline]

sets/returns z componet of the magnetic field.

Definition at line 289 of file Helix.h.

References m_alpha, m_bField, and updateCache().

00289                        {
00290     m_bField = a;
00291     m_alpha = 10000. / 2.99792458 / m_bField;
00292     updateCache();
00293     return m_bField;
00294 }

const HepPoint3D & Helix::center ( void   )  const [inline]

returns position of helix center(z = 0.);

Definition at line 203 of file Helix.h.

References m_center.

Referenced by TTrack::approach(), TTrack::approach2D(), TBuilderCurl::buildStereoMC(), calVirtualCircle(), TTrack::center(), TCurlFinder::distance(), TTrack::dump(), TTrack::dxda2D(), Emc_helix::Emc_Get(), EsTimeAlg::execute(), TConformalFinder0::findCloseClusters(), TFastFinder::findCloseHits(), TCurlFinder::findCloseHits(), TConformalFinder0::findCloseHits(), TTrack::fit2D(), TTrack::HelCyl(), TTrack::impact(), TRunge::intersect_cylinder(), TRunge::intersect_yz_plane(), TRunge::intersect_zx_plane(), TCurlFinder::mask3DTrack(), TTrackManager::merge(), TCurlFinder::merge3DTrack(), TCurlFinder::plotTrack(), TCurlFinder::salvage3DTrack(), TTrack::stereoHitForCurl(), TTrack::szPosition(), TofFz_helix::TofFz_Get(), and TCurlFinder::trace3DTrack().

00203                         {
00204     return m_center;
00205 }

double Helix::cosPhi0 ( void   )  const [inline]

Definition at line 310 of file Helix.h.

References m_cp.

00310                          {
00311     return m_cp;
00312 }

double Helix::cosTheta ( void   )  const [inline]

Definition at line 126 of file Helix.h.

References m_a.

Referenced by MdcHoughFinder::GetMcInfo(), and HoughValidUpdate::GetMcInfo().

00126 { return m_a[4]/sqrt(1.+ m_a[4]*m_a[4]); }

double Helix::curv ( void   )  const [inline]

Definition at line 257 of file Helix.h.

References m_r.

Referenced by TBuilder0::buildStereo(), TBuilder0::buildStereo0(), TBuilderCurl::buildStereoMC(), TFastFinder::findCloseHits(), TCurlFinder::findCloseHits(), TConformalFinder0::findCloseHits(), TTrack::stereoHitForCurl(), and TTrack::szPosition().

00257                       {
00258     return m_r;
00259 }

HepMatrix Helix::del4MDelA ( double  phi,
double  mass 
) const

Definition at line 597 of file Helix.cxx.

References cos(), DBL_MAX, m_ac, phi0(), and sin().

Referenced by momentum(), and KalmanFit::Helix::momentum().

00597                                               {
00598     //
00599     //   Calculate Jacobian (@4m/@a)
00600     //   Vector a  is helix parameters and phi is internal parameter.
00601     //   Vector 4m is 4 momentum.
00602     //
00603 
00604     HepMatrix d4MDA(4,5,0);
00605 
00606     double phi0 = m_ac[1];
00607     double cpa  = m_ac[2];
00608     double tnl  = m_ac[4];
00609 
00610     double cosf0phi = cos(phi0+phi);
00611     double sinf0phi = sin(phi0+phi);
00612 
00613     double rho;
00614     if(cpa != 0.)rho = 1./cpa;
00615     else rho = (DBL_MAX);
00616 
00617     double charge = 1.;
00618     if(cpa < 0.)charge = -1.;
00619 
00620     double E = sqrt(rho*rho*(1.+tnl*tnl)+mass*mass);
00621 
00622     d4MDA[0][1] = -fabs(rho)*cosf0phi;
00623     d4MDA[0][2] = charge*rho*rho*sinf0phi;
00624  
00625     d4MDA[1][1] = -fabs(rho)*sinf0phi;
00626     d4MDA[1][2] = -charge*rho*rho*cosf0phi;
00627 
00628     d4MDA[2][2] = -charge*rho*rho*tnl;
00629     d4MDA[2][4] = fabs(rho);
00630 
00631     if (cpa != 0.0 && E != 0.0) {
00632       d4MDA[3][2] = (-1.-tnl*tnl)/(cpa*cpa*cpa*E);
00633       d4MDA[3][4] = tnl/(cpa*cpa*E);
00634     } else {
00635       d4MDA[3][2] = (DBL_MAX);
00636       d4MDA[3][4] = (DBL_MAX);
00637     }
00638     return d4MDA;
00639 }

HepMatrix Helix::del4MXDelA ( double  phi,
double  mass 
) const

Definition at line 643 of file Helix.cxx.

References cos(), DBL_MAX, dr(), dz(), m_ac, m_cp, m_r, m_sp, phi0(), and sin().

Referenced by momentum(), and KalmanFit::Helix::momentum().

00643                                                {
00644     //
00645     //   Calculate Jacobian (@4mx/@a)
00646     //   Vector a  is helix parameters and phi is internal parameter.
00647     //   Vector 4xm is 4 momentum and position.
00648     //
00649 
00650     HepMatrix d4MXDA(7,5,0);
00651 
00652     const double & dr      = m_ac[0];
00653     const double & phi0    = m_ac[1];
00654     const double & cpa     = m_ac[2];
00655     const double & dz      = m_ac[3];
00656     const double & tnl     = m_ac[4];
00657 
00658     double cosf0phi = cos(phi0+phi);
00659     double sinf0phi = sin(phi0+phi);
00660 
00661     double rho;
00662     if(cpa != 0.)rho = 1./cpa;
00663     else rho = (DBL_MAX);
00664 
00665     double charge = 1.;
00666     if(cpa < 0.)charge = -1.;
00667 
00668     double E = sqrt(rho * rho * (1. + tnl * tnl) + mass * mass);
00669 
00670     d4MXDA[0][1] = - fabs(rho) * cosf0phi;
00671     d4MXDA[0][2] = charge * rho * rho * sinf0phi;
00672  
00673     d4MXDA[1][1] = - fabs(rho) * sinf0phi;
00674     d4MXDA[1][2] = - charge * rho * rho * cosf0phi;
00675 
00676     d4MXDA[2][2] = - charge * rho * rho * tnl;
00677     d4MXDA[2][4] = fabs(rho);
00678 
00679     if (cpa != 0.0 && E != 0.0) {
00680       d4MXDA[3][2] = (- 1. - tnl * tnl) / (cpa * cpa * cpa * E);
00681       d4MXDA[3][4] = tnl / (cpa * cpa * E);
00682     } else {
00683       d4MXDA[3][2] = (DBL_MAX);
00684       d4MXDA[3][4] = (DBL_MAX);
00685     }
00686     
00687     d4MXDA[4][0] = m_cp;
00688     d4MXDA[4][1] = - dr * m_sp + m_r * (- m_sp + sinf0phi);
00689     if (cpa != 0.0) {
00690       d4MXDA[4][2] = - (m_r / cpa) * (m_cp - cosf0phi);
00691     } else {
00692       d4MXDA[4][2] = (DBL_MAX);
00693     }
00694      
00695     d4MXDA[5][0] = m_sp;
00696     d4MXDA[5][1] = dr * m_cp + m_r * (m_cp - cosf0phi);
00697     if (cpa != 0.0) {
00698       d4MXDA[5][2] = - (m_r / cpa) * (m_sp - sinf0phi);
00699 
00700       d4MXDA[6][2] = (m_r / cpa) * tnl * phi;
00701     } else {
00702       d4MXDA[5][2] = (DBL_MAX);
00703 
00704       d4MXDA[6][2] = (DBL_MAX);
00705     }
00706 
00707     d4MXDA[6][3] = 1.;
00708     d4MXDA[6][4] = - m_r * phi;
00709 
00710     return d4MXDA;
00711 }

HepMatrix Helix::delApDelA ( const HepVector &  ap  )  const

Definition at line 438 of file Helix.cxx.

References cos(), DBL_MAX, dr(), dz(), m_ac, M_PI, M_PI2, M_PI8, m_r, phi0(), and sin().

Referenced by pivot(), and KalmanFit::Helix::pivot().

00438                                            {
00439     //
00440     //   Calculate Jacobian (@ap/@a)
00441     //   Vector ap is new helix parameters and a is old helix parameters. 
00442     //
00443 
00444     HepMatrix dApDA(5,5,0);
00445 
00446     const double & dr    = m_ac[0];
00447     const double & phi0  = m_ac[1];
00448     const double & cpa   = m_ac[2];
00449     const double & dz    = m_ac[3];
00450     const double & tnl   = m_ac[4];
00451 
00452     double drp   = ap[0];
00453     double phi0p = ap[1];
00454     double cpap  = ap[2];
00455     double dzp   = ap[3];
00456     double tnlp  = ap[4];
00457 
00458     double rdr   = m_r + dr;
00459     double rdrpr;
00460     if ((m_r + drp) != 0.0) {
00461       rdrpr = 1. / (m_r + drp);
00462     } else {
00463       rdrpr = (DBL_MAX);
00464     }
00465     // double csfd  = cos(phi0)*cos(phi0p) + sin(phi0)*sin(phi0p); 
00466     // double snfd  = cos(phi0)*sin(phi0p) - sin(phi0)*cos(phi0p); 
00467     double csfd  = cos(phi0p - phi0);
00468     double snfd  = sin(phi0p - phi0);
00469     double phid  = fmod(phi0p - phi0 + M_PI8, M_PI2);
00470     if (phid > M_PI) phid = phid - M_PI2;
00471 
00472     dApDA[0][0]  =  csfd;
00473     dApDA[0][1]  =  rdr*snfd;
00474     if(cpa!=0.0) {
00475       dApDA[0][2]  =  (m_r/cpa)*( 1.0 - csfd );
00476     } else {
00477       dApDA[0][2]  = (DBL_MAX);
00478     }
00479 
00480     dApDA[1][0]  = - rdrpr*snfd;
00481     dApDA[1][1]  = rdr*rdrpr*csfd;
00482     if(cpa!=0.0) {
00483       dApDA[1][2]  = (m_r/cpa)*rdrpr*snfd;
00484     } else {
00485       dApDA[1][2]  = (DBL_MAX);
00486     }
00487     
00488     dApDA[2][2]  = 1.0;
00489 
00490     dApDA[3][0]  = m_r*rdrpr*tnl*snfd;
00491     dApDA[3][1]  = m_r*tnl*(1.0 - rdr*rdrpr*csfd);
00492     if(cpa!=0.0) {
00493       dApDA[3][2]  = (m_r/cpa)*tnl*(phid - m_r*rdrpr*snfd);
00494     } else {
00495       dApDA[3][2]  = (DBL_MAX);
00496     }
00497     dApDA[3][3]  = 1.0;
00498     dApDA[3][4]  = - m_r*phid;
00499 
00500     dApDA[4][4] = 1.0;
00501  
00502     return dApDA;
00503 }

HepMatrix Helix::delMDelA ( double  phi  )  const

Definition at line 560 of file Helix.cxx.

References cos(), DBL_MAX, m_ac, phi0(), and sin().

Referenced by momentum(), and KalmanFit::Helix::momentum().

00560                                 {
00561     //
00562     //   Calculate Jacobian (@m/@a)
00563     //   Vector a is helix parameters and phi is internal parameter.
00564     //   Vector m is momentum.
00565     //
00566 
00567     HepMatrix dMDA(3,5,0);
00568 
00569     const double & phi0 = m_ac[1];
00570     const double & cpa  = m_ac[2];
00571     const double & tnl  = m_ac[4];
00572 
00573     double cosf0phi = cos(phi0+phi);
00574     double sinf0phi = sin(phi0+phi);
00575 
00576     double rho;
00577     if(cpa != 0.)rho = 1./cpa;
00578     else rho = (DBL_MAX);
00579 
00580     double charge = 1.;
00581     if(cpa < 0.)charge = -1.;
00582 
00583     dMDA[0][1] = -fabs(rho)*cosf0phi;
00584     dMDA[0][2] = charge*rho*rho*sinf0phi;
00585  
00586     dMDA[1][1] = -fabs(rho)*sinf0phi;
00587     dMDA[1][2] = -charge*rho*rho*cosf0phi;
00588 
00589     dMDA[2][2] = -charge*rho*rho*tnl;
00590     dMDA[2][4] = fabs(rho);
00591 
00592     return dMDA;
00593 }

HepMatrix Helix::delXDelA ( double  phi  )  const

Definition at line 506 of file Helix.cxx.

References cos(), DBL_MAX, dr(), dz(), m_ac, m_cp, m_r, m_sp, phi0(), and sin().

Referenced by x(), and KalmanFit::Helix::x().

00506                                 {
00507     //
00508     //   Calculate Jacobian (@x/@a)
00509     //   Vector a is helix parameters and phi is internal parameter
00510     //   which specifys the point to be calculated for Ex(phi).
00511     //
00512 
00513     HepMatrix dXDA(3,5,0);
00514 
00515     const double & dr      = m_ac[0];
00516     const double & phi0    = m_ac[1];
00517     const double & cpa     = m_ac[2];
00518     const double & dz      = m_ac[3];
00519     const double & tnl     = m_ac[4];
00520 
00521     double cosf0phi = cos(phi0 + phi);
00522     double sinf0phi = sin(phi0 + phi);
00523 
00524     dXDA[0][0]     = m_cp;
00525     dXDA[0][1]     = - dr * m_sp + m_r * (- m_sp + sinf0phi); 
00526     if(cpa!=0.0) {
00527       dXDA[0][2]     = - (m_r / cpa) * (m_cp - cosf0phi);
00528     } else {
00529       dXDA[0][2] = (DBL_MAX);
00530     }
00531     // dXDA[0][3]     = 0.0;
00532     // dXDA[0][4]     = 0.0;
00533  
00534     dXDA[1][0]     = m_sp;
00535     dXDA[1][1]     = dr * m_cp + m_r * (m_cp - cosf0phi);
00536     if(cpa!=0.0) {
00537       dXDA[1][2]     = - (m_r / cpa) * (m_sp - sinf0phi);
00538     } else {
00539       dXDA[1][2] = (DBL_MAX);
00540     }
00541     // dXDA[1][3]     = 0.0;
00542     // dXDA[1][4]     = 0.0;
00543 
00544     // dXDA[2][0]     = 0.0;
00545     // dXDA[2][1]     = 0.0;
00546     if(cpa!=0.0) {
00547       dXDA[2][2]     = (m_r / cpa) * tnl * phi;
00548     } else {
00549       dXDA[2][2] = (DBL_MAX);
00550     }
00551     dXDA[2][3]     = 1.0;
00552     dXDA[2][4]     = - m_r * phi;
00553 
00554     return dXDA;
00555 }

Hep3Vector Helix::direction ( double  dPhi = 0.  )  const [inline]

returns direction vector after rotating angle dPhi in phi direction.

Definition at line 221 of file Helix.h.

References momentum().

Referenced by Emc_helix::Emc_Get(), THelixFitter::main(), and TofFz_helix::TofFz_Get().

00221                                  {
00222     return momentum(phi).unit();
00223 }

double Helix::dr ( void   )  const [inline]

returns an element of parameters.

Definition at line 227 of file Helix.h.

References m_ac.

Referenced by del4MXDelA(), KalmanFit::Helix::del4MXDelA(), delApDelA(), KalmanFit::Helix::delApDelA(), delXDelA(), KalmanFit::Helix::delXDelA(), MdcHoughFinder::GetMcInfo(), HoughValidUpdate::GetMcInfo(), FTFinder::makeMdst(), pivot(), and KalmanFit::Helix::pivot().

00227                     {
00228     return m_ac[0];
00229 }

double Helix::dz ( void   )  const [inline]

Definition at line 245 of file Helix.h.

References m_ac.

Referenced by del4MXDelA(), KalmanFit::Helix::del4MXDelA(), delApDelA(), KalmanFit::Helix::delApDelA(), delXDelA(), KalmanFit::Helix::delXDelA(), MdcHoughFinder::GetMcInfo(), HoughValidUpdate::GetMcInfo(), TCurlFinder::makeCurlTracks(), FTFinder::makeMdst(), TTrackManager::maskCurl(), TTrackManager::merge(), TCurlFinder::merge3DTrack(), pivot(), KalmanFit::Helix::pivot(), and TRunge::SetFlightLength().

00245                     {
00246     return m_ac[3];
00247 }

const HepSymMatrix & Helix::Ea ( const HepSymMatrix &  newdA  )  [inline]

sets helix paramters and error matrix.

Definition at line 283 of file Helix.h.

References m_Ea.

00283                                 {
00284     return m_Ea = i;
00285 }

const HepSymMatrix & Helix::Ea ( void   )  const [inline]

returns error matrix.

Definition at line 269 of file Helix.h.

References m_Ea.

Referenced by TTrackMC::compare(), TTrack::fit2D(), HelixHasNan(), FTFinder::makeTds(), TRunge::pivot(), and PositiveDefinite().

00269                     {
00270     return m_Ea;
00271 }

void Helix::ignoreErrorMatrix ( void   ) 

unsets error matrix. Error calculations will be ignored after this function call until an error matrix be set again. 0 matrix will be return as a return value for error matrix when you call functions which returns an error matrix.

Definition at line 714 of file Helix.cxx.

References m_Ea, and m_matrixValid.

Referenced by EsTimeAlg::execute().

00714                              {
00715     m_matrixValid = false;
00716     m_Ea *= 0.;
00717 }

double Helix::kappa ( void   )  const [inline]

Definition at line 239 of file Helix.h.

References m_ac.

Referenced by TTrack::approach(), TTrackMC::compare(), FTFinder::makeMdst(), pivot(), KalmanFit::Helix::pivot(), and TTrack::TTrack().

00239                        {
00240     return m_ac[2];
00241 }

HepLorentzVector Helix::momentum ( double  dPhi,
double  mass,
HepPoint3D x,
HepSymMatrix &  Emx 
) const

returns 4momentum vector after rotating angle dPhi in phi direction.

Definition at line 275 of file Helix.cxx.

References cos(), del4MXDelA(), m_ac, m_cp, m_Ea, m_matrixValid, m_pivot, m_pt, m_r, m_sp, pt(), and sin().

00278                                           {
00279   // 
00280   // Calculate momentum.
00281   //
00282   // Pt = | 1/kappa | (GeV/c)
00283   //
00284   // Px = -Pt * sin(phi0 + phi) 
00285   // Py =  Pt * cos(phi0 + phi)
00286   // Pz =  Pt * tan(lambda)
00287   //
00288   // E  = sqrt( 1/kappa/kappa * (1+tan(lambda)*tan(lambda)) + mass*mass )
00289   
00290   double pt = fabs(m_pt);
00291   double px = - pt * sin(m_ac[1] + phi);
00292   double py =   pt * cos(m_ac[1] + phi);
00293   double pz =   pt * m_ac[4];
00294   double E  = sqrt(pt * pt * (1. + m_ac[4] * m_ac[4]) + mass * mass);
00295 
00296   x.setX(m_pivot.x() + m_ac[0] * m_cp + m_r * (m_cp - cos(m_ac[1] + phi)));
00297   x.setY(m_pivot.y() + m_ac[0] * m_sp + m_r * (m_sp - sin(m_ac[1] + phi)));
00298   x.setZ(m_pivot.z() + m_ac[3] - m_r * m_ac[4] * phi);
00299 
00300   if (m_matrixValid) Emx = m_Ea.similarity(del4MXDelA(phi,mass));
00301   else               Emx = m_Ea;
00302   
00303   return HepLorentzVector(px, py, pz, E);
00304 }

HepLorentzVector Helix::momentum ( double  dPhi,
double  mass,
HepSymMatrix &  Em 
) const

returns 4momentum vector after rotating angle dPhi in phi direction.

Definition at line 250 of file Helix.cxx.

References cos(), del4MDelA(), m_ac, m_Ea, m_matrixValid, m_pt, pt(), and sin().

00250                                                                 {
00251     // 
00252     // Calculate momentum.
00253     //
00254     // Pt = | 1/kappa | (GeV/c)
00255     //
00256     // Px = -Pt * sin(phi0 + phi) 
00257     // Py =  Pt * cos(phi0 + phi)
00258     // Pz =  Pt * tan(lambda)
00259     //
00260     // E  = sqrt( 1/kappa/kappa * (1+tan(lambda)*tan(lambda)) + mass*mass )
00261 
00262     double pt = fabs(m_pt);
00263     double px = - pt * sin(m_ac[1] + phi);
00264     double py =   pt * cos(m_ac[1] + phi);
00265     double pz =   pt * m_ac[4];
00266     double E  =   sqrt(pt*pt*(1.+m_ac[4]*m_ac[4])+mass*mass);
00267 
00268     if (m_matrixValid) Em = m_Ea.similarity(del4MDelA(phi,mass));
00269     else               Em = m_Ea;
00270 
00271     return HepLorentzVector(px, py, pz, E);
00272 }

HepLorentzVector Helix::momentum ( double  dPhi,
double  mass 
) const

returns 4momentum vector after rotating angle dPhi in phi direction.

Definition at line 227 of file Helix.cxx.

References cos(), m_ac, m_pt, pt(), and sin().

00227                                              {
00228     // 
00229     // Calculate momentum.
00230     //
00231     // Pt = | 1/kappa | (GeV/c)
00232     //
00233     // Px = -Pt * sin(phi0 + phi) 
00234     // Py =  Pt * cos(phi0 + phi)
00235     // Pz =  Pt * tan(lambda)
00236     //
00237     // E  = sqrt( 1/kappa/kappa * (1+tan(lambda)*tan(lambda)) + mass*mass )
00238 
00239     double pt = fabs(m_pt);
00240     double px = - pt * sin(m_ac[1] + phi);
00241     double py =   pt * cos(m_ac[1] + phi);
00242     double pz =   pt * m_ac[4];
00243     double E  =   sqrt(pt*pt*(1.+m_ac[4]*m_ac[4])+mass*mass);
00244 
00245     return HepLorentzVector(px, py, pz, E);
00246 }

Hep3Vector Helix::momentum ( double  dPhi,
HepSymMatrix &  Em 
) const

returns momentum vector after rotating angle dPhi in phi direction.

Definition at line 204 of file Helix.cxx.

References cos(), delMDelA(), m_ac, m_Ea, m_matrixValid, m_pt, pt(), and sin().

00204                                                    {
00205     // 
00206     // Calculate momentum.
00207     //
00208     // Pt = | 1/kappa | (GeV/c)
00209     //
00210     // Px = -Pt * sin(phi0 + phi) 
00211     // Py =  Pt * cos(phi0 + phi)
00212     // Pz =  Pt * tan(lambda)
00213     //
00214 
00215     double pt = fabs(m_pt);
00216     double px = - pt * sin(m_ac[1] + phi);
00217     double py =   pt * cos(m_ac[1] + phi);
00218     double pz =   pt * m_ac[4];
00219 
00220     if (m_matrixValid) Em = m_Ea.similarity(delMDelA(phi));
00221     else               Em = m_Ea;
00222 
00223     return Hep3Vector(px, py, pz);
00224 }

Hep3Vector Helix::momentum ( double  dPhi = 0.  )  const

returns momentum vector after rotating angle dPhi in phi direction.

Definition at line 184 of file Helix.cxx.

References cos(), m_ac, m_pt, pt(), and sin().

Referenced by TTrackMC::compare(), direction(), KalmanFit::Helix::direction(), EsTimeAlg::execute(), MdcHoughFinder::GetMcInfo(), HoughValidUpdate::GetMcInfo(), FTFinder::makeMdst(), and TTrack::p().

00184                                 {
00185     // 
00186     // Calculate momentum.
00187     //
00188     // Pt = | 1/kappa | (GeV/c)
00189     //
00190     // Px = -Pt * sin(phi0 + phi) 
00191     // Py =  Pt * cos(phi0 + phi)
00192     // Pz =  Pt * tan(lambda)
00193     //
00194 
00195     double pt = fabs(m_pt);
00196     double px = - pt * sin(m_ac[1] + phi);
00197     double py =   pt * cos(m_ac[1] + phi);
00198     double pz =   pt * m_ac[4];
00199 
00200     return Hep3Vector(px, py, pz);
00201 }

Helix & Helix::operator= ( const Helix  ) 

Copy operator.

Definition at line 380 of file Helix.cxx.

References genRecEmupikp::i, m_a, m_ac, m_alpha, m_bField, m_center, m_cp, m_Ea, m_matrixValid, m_pivot, m_pt, m_r, and m_sp.

00380                                   {
00381     if (this == & i) return * this;
00382 
00383     m_bField = i.m_bField;
00384     m_alpha = i.m_alpha;
00385     m_pivot = i.m_pivot;
00386     m_a = i.m_a;
00387     m_Ea = i.m_Ea;
00388     m_matrixValid = i.m_matrixValid;
00389 
00390     m_center = i.m_center;
00391     m_cp = i.m_cp;
00392     m_sp = i.m_sp;
00393     m_pt = i.m_pt;
00394     m_r  = i.m_r;
00395     m_ac[0] = i.m_ac[0];
00396     m_ac[1] = i.m_ac[1];
00397     m_ac[2] = i.m_ac[2];
00398     m_ac[3] = i.m_ac[3];
00399     m_ac[4] = i.m_ac[4];
00400 
00401     return * this;
00402 }

double Helix::phi0 ( void   )  const [inline]

Definition at line 233 of file Helix.h.

References m_ac.

Referenced by TTrack::approach(), del4MDelA(), KalmanFit::Helix::del4MDelA(), del4MXDelA(), KalmanFit::Helix::del4MXDelA(), delApDelA(), KalmanFit::Helix::delApDelA(), delMDelA(), KalmanFit::Helix::delMDelA(), delXDelA(), KalmanFit::Helix::delXDelA(), MdcHoughFinder::GetMcInfo(), HoughValidUpdate::GetMcInfo(), TRunge::intersect_cylinder(), FTFinder::makeMdst(), pivot(), and KalmanFit::Helix::pivot().

00233                       {
00234     return m_ac[1];
00235 }

const HepPoint3D & Helix::pivot ( const HepPoint3D newPivot  ) 

sets pivot position.

Definition at line 308 of file Helix.cxx.

References cos(), delApDelA(), dr(), dz(), kappa(), m_a, m_ac, m_Ea, m_matrixValid, M_PI, M_PI2, M_PI4, M_PI8, m_pivot, m_r, phi0(), tanl(), and updateCache().

00308                                         {
00309     const double & dr    = m_ac[0];
00310     const double & phi0  = m_ac[1];
00311     const double & kappa = m_ac[2];
00312     const double & dz    = m_ac[3];
00313     const double & tanl  = m_ac[4];
00314 
00315     double rdr = dr + m_r;
00316     double phi = fmod(phi0 + M_PI4, M_PI2);
00317     double csf0 = cos(phi);
00318     double snf0 = (1. - csf0) * (1. + csf0);
00319     snf0 = sqrt((snf0 > 0.) ? snf0 : 0.);
00320     if(phi > M_PI) snf0 = - snf0;
00321 
00322     double xc = m_pivot.x() + rdr * csf0;
00323     double yc = m_pivot.y() + rdr * snf0;
00324     double csf, snf;
00325     if(m_r != 0.0) {
00326       csf = (xc - newPivot.x()) / m_r;
00327       snf = (yc - newPivot.y()) / m_r;
00328       double anrm = sqrt(csf * csf + snf * snf);
00329       if(anrm != 0.0) {
00330         csf /= anrm;
00331         snf /= anrm;
00332         phi = atan2(snf, csf);
00333       } else {
00334         csf = 1.0;
00335         snf = 0.0;
00336         phi = 0.0;
00337       }
00338     } else {
00339       csf = 1.0;
00340       snf = 0.0;
00341       phi = 0.0;
00342     }
00343     double phid = fmod(phi - phi0 + M_PI8, M_PI2);
00344     if(phid > M_PI) phid = phid - M_PI2;
00345     double drp = (m_pivot.x() + dr * csf0 + m_r * (csf0 - csf) - newPivot.x())
00346         * csf
00347         + (m_pivot.y() + dr * snf0 + m_r * (snf0 - snf) - newPivot.y()) * snf;
00348     double dzp = m_pivot.z() + dz - m_r * tanl * phid - newPivot.z();
00349 
00350     HepVector ap(5);
00351     ap[0] = drp;
00352     ap[1] = fmod(phi + M_PI4, M_PI2);
00353     ap[2] = kappa;
00354     ap[3] = dzp;
00355     ap[4] = tanl;
00356 
00357     //    if (m_matrixValid) m_Ea.assign(delApDelA(ap) * m_Ea * delApDelA(ap).T());
00358     if (m_matrixValid) m_Ea = m_Ea.similarity(delApDelA(ap));
00359 
00360     m_a = ap;
00361     m_pivot = newPivot;
00362 
00363     //...Are these needed?...iw...
00364     updateCache();
00365     return m_pivot;
00366 }

const HepPoint3D & Helix::pivot ( void   )  const [inline]

returns pivot position.

Definition at line 209 of file Helix.h.

References m_pivot.

Referenced by TTrackMC::compare(), TPerfectFinder::doit(), TTrack::dump(), Emc_helix::Emc_Get(), EsTimeAlg::execute(), MdcHoughFinder::GetMcInfo(), HoughValidUpdate::GetMcInfo(), TTrack::HelCyl(), FTFinder::makeTds(), TTrackManager::maskCurl(), TTrackManager::merge(), TTrack::movePivot(), TRunge::pivot(), TRunge::SetFlightLength(), TofFz_helix::TofFz_Get(), and TrackKinematics().

00209                        {
00210     return m_pivot;
00211 }

double Helix::pt ( void   )  const [inline]

Definition at line 125 of file Helix.h.

References m_pt.

Referenced by MdcHoughFinder::GetMcInfo(), HoughValidUpdate::GetMcInfo(), momentum(), and KalmanFit::Helix::momentum().

00125 { return m_pt; }

double Helix::radius ( void   )  const [inline]

returns radious of helix.

Definition at line 215 of file Helix.h.

References m_r.

Referenced by TCurlFinder::distance(), Emc_helix::Emc_Get(), EsTimeAlg::execute(), TConformalFinder0::findCloseClusters(), TTrack::fit2D(), HoughValidUpdate::GetMcInfo(), TTrack::HelCyl(), TTrack::impact(), RkFitCylinder::intersect(), TRunge::intersect_cylinder(), TRunge::intersect_xy_plane(), TRunge::intersect_yz_plane(), TRunge::intersect_zx_plane(), TCurlFinder::mask3DTrack(), TTrackManager::merge(), TCurlFinder::merge3DTrack(), TCurlFinder::plotTrack(), TTrack::radius(), TofFz_helix::TofFz_Get(), and TCurlFinder::trace3DTrack().

00215                         {
00216     return m_r;
00217 }

void Helix::set ( const HepPoint3D pivot,
const HepVector &  a,
const HepSymMatrix &  Ea 
)

sets helix pivot position, parameters, and error matrix.

Definition at line 369 of file Helix.cxx.

References m_a, m_Ea, m_matrixValid, m_pivot, and updateCache().

00371                                     {
00372     m_pivot = pivot;
00373     m_a = a;
00374     m_Ea = Ea;
00375     m_matrixValid = true;
00376     updateCache();
00377 }

double Helix::sinPhi0 ( void   )  const [inline]

Definition at line 304 of file Helix.h.

References m_sp.

00304                          {
00305     return m_sp;
00306 }

double Helix::tanl ( void   )  const [inline]

Definition at line 251 of file Helix.h.

References m_ac.

Referenced by TTrack::approach(), TTrackMC::compare(), MdcHoughFinder::GetMcInfo(), HoughValidUpdate::GetMcInfo(), RkFitCylinder::intersect(), TRunge::intersect_xy_plane(), FTFinder::makeMdst(), FTFinder::makeTds(), pivot(), and KalmanFit::Helix::pivot().

00251                       {
00252     return m_ac[4];
00253 }

void Helix::updateCache ( void   )  [private]

Definition at line 405 of file Helix.cxx.

References cos(), DBL_MAX, m_a, m_ac, m_alpha, m_center, m_cp, m_pivot, m_pt, m_r, m_sp, sin(), and x().

Referenced by a(), KalmanFit::Helix::a(), bFieldZ(), KalmanFit::Helix::bFieldZ(), Helix(), pivot(), KalmanFit::Helix::pivot(), set(), and KalmanFit::Helix::set().

00405                        {
00406     //
00407     //   Calculate Helix center( xc, yc ).
00408     //   
00409     //   xc = x0 + (dr + (alpha / kappa)) * cos(phi0)  (cm)
00410     //   yc = y0 + (dr + (alpha / kappa)) * sin(phi0)  (cm)
00411     //
00412 
00413     m_ac[0] = m_a[0];
00414     m_ac[1] = m_a[1];
00415     m_ac[2] = m_a[2];
00416     m_ac[3] = m_a[3];
00417     m_ac[4] = m_a[4];
00418 
00419     m_cp = cos(m_ac[1]);
00420     m_sp = sin(m_ac[1]);
00421     if (m_ac[2] != 0.0) {
00422         m_pt = 1. / m_ac[2];
00423         m_r = m_alpha / m_ac[2];
00424     }
00425     else {
00426         m_pt = (DBL_MAX);
00427         m_r = (DBL_MAX);
00428     }
00429 
00430     double x = m_pivot.x() + (m_ac[0] + m_r) * m_cp;
00431     double y = m_pivot.y() + (m_ac[0] + m_r) * m_sp;
00432     m_center.setX(x);
00433     m_center.setY(y);
00434     m_center.setZ(0.);
00435 }

HepPoint3D Helix::x ( double  dPhi,
HepSymMatrix &  Ex 
) const

returns position and convariance matrix(Ex) after rotation.

Definition at line 163 of file Helix.cxx.

References cos(), delXDelA(), m_ac, m_cp, m_Ea, m_matrixValid, m_pivot, m_r, m_sp, sin(), and x().

00163                                             {
00164     double x = m_pivot.x() + m_ac[0] * m_cp + m_r * (m_cp - cos(m_ac[1] +phi));
00165     double y = m_pivot.y() + m_ac[0] * m_sp + m_r * (m_sp - sin(m_ac[1] +phi));
00166     double z = m_pivot.z() + m_ac[3] - m_r * m_ac[4] * phi;
00167 
00168     //
00169     //   Calculate position error matrix.
00170     //   Ex(phi) = (@x/@a)(Ea)(@x/@a)^T, phi is deflection angle to specify the
00171     //   point to be calcualted.
00172     //
00173     // HepMatrix dXDA(3, 5, 0);
00174     // dXDA = delXDelA(phi);
00175     // Ex.assign(dXDA * m_Ea * dXDA.T());
00176 
00177     if (m_matrixValid) Ex = m_Ea.similarity(delXDelA(phi));
00178     else               Ex = m_Ea;
00179 
00180     return HepPoint3D(x, y, z);
00181 }

double * Helix::x ( double  dPhi,
double  p[3] 
) const

Definition at line 146 of file Helix.cxx.

References cos(), m_ac, m_cp, m_pivot, m_r, m_sp, and sin().

00146                                       {
00147     //
00148     // Calculate position (x,y,z) along helix.
00149     //   
00150     // x = x0 + dr * cos(phi0) + (alpha / kappa) * (cos(phi0) - cos(phi0+phi))
00151     // y = y0 + dr * sin(phi0) + (alpha / kappa) * (sin(phi0) - sin(phi0+phi))
00152     // z = z0 + dz             - (alpha / kappa) * tan(lambda) * phi
00153     //
00154 
00155     p[0] = m_pivot.x() + m_ac[0] * m_cp + m_r * (m_cp - cos(m_ac[1] + phi));
00156     p[1] = m_pivot.y() + m_ac[0] * m_sp + m_r * (m_sp - sin(m_ac[1] + phi));
00157     p[2] = m_pivot.z() + m_ac[3] - m_r * m_ac[4] * phi;
00158 
00159     return p;
00160 }

HepPoint3D Helix::x ( double  dPhi = 0.  )  const

returns position after rotating angle dPhi in phi direction.

Definition at line 129 of file Helix.cxx.

References cos(), m_ac, m_cp, m_pivot, m_r, m_sp, and sin().

Referenced by TTrack::approach(), TTrack::approach2D(), TCurlFinder::distance(), TTrack::dxda(), EsTimeAlg::execute(), TTrack::fit2D(), TTrack::HelCyl(), RkFitCylinder::intersect(), FTFinder::makeMdst(), TRunge::SetFlightLength(), TofFz_helix::TofFz_Get(), updateCache(), KalmanFit::Helix::updateCache(), x(), and KalmanFit::Helix::x().

00129                          {
00130     //
00131     // Calculate position (x,y,z) along helix.
00132     //   
00133     // x = x0 + dr * cos(phi0) + (alpha / kappa) * (cos(phi0) - cos(phi0+phi))
00134     // y = y0 + dr * sin(phi0) + (alpha / kappa) * (sin(phi0) - sin(phi0+phi))
00135     // z = z0 + dz             - (alpha / kappa) * tan(lambda) * phi
00136     //
00137 
00138     double x = m_pivot.x() + m_ac[0] * m_cp + m_r * (m_cp - cos(m_ac[1] +phi));
00139     double y = m_pivot.y() + m_ac[0] * m_sp + m_r * (m_sp - sin(m_ac[1] +phi));
00140     double z = m_pivot.z() + m_ac[3] - m_r * m_ac[4] * phi;
00141 
00142     return HepPoint3D(x, y, z);
00143 }


Member Data Documentation

const double Helix::ConstantAlpha = 333.564095 [static]

Constant alpha for uniform field.

Definition at line 171 of file Helix.h.

HepVector Helix::m_a [private]

Definition at line 177 of file Helix.h.

Referenced by a(), KalmanFit::Helix::a(), cosTheta(), Helix(), operator=(), KalmanFit::Helix::operator=(), pivot(), KalmanFit::Helix::pivot(), set(), KalmanFit::Helix::set(), updateCache(), and KalmanFit::Helix::updateCache().

double Helix::m_ac[5] [private]

Definition at line 187 of file Helix.h.

Referenced by del4MDelA(), KalmanFit::Helix::del4MDelA(), del4MXDelA(), KalmanFit::Helix::del4MXDelA(), delApDelA(), KalmanFit::Helix::delApDelA(), delMDelA(), KalmanFit::Helix::delMDelA(), delXDelA(), KalmanFit::Helix::delXDelA(), dr(), KalmanFit::Helix::dr(), dz(), KalmanFit::Helix::dz(), kappa(), KalmanFit::Helix::kappa(), momentum(), KalmanFit::Helix::momentum(), operator=(), KalmanFit::Helix::operator=(), phi0(), KalmanFit::Helix::phi0(), pivot(), KalmanFit::Helix::pivot(), tanl(), KalmanFit::Helix::tanl(), updateCache(), KalmanFit::Helix::updateCache(), x(), and KalmanFit::Helix::x().

double Helix::m_alpha [protected]

Definition at line 164 of file Helix.h.

Referenced by KalmanFit::Helix::alpha(), bFieldZ(), KalmanFit::Helix::bFieldZ(), Helix(), operator=(), KalmanFit::Helix::operator=(), updateCache(), and KalmanFit::Helix::updateCache().

double Helix::m_bField [protected]

Definition at line 163 of file Helix.h.

Referenced by bFieldZ(), KalmanFit::Helix::bFieldZ(), Helix(), operator=(), and KalmanFit::Helix::operator=().

HepPoint3D Helix::m_center [private]

Definition at line 182 of file Helix.h.

Referenced by center(), KalmanFit::Helix::center(), operator=(), KalmanFit::Helix::operator=(), updateCache(), and KalmanFit::Helix::updateCache().

double Helix::m_cp [private]

Definition at line 183 of file Helix.h.

Referenced by cosPhi0(), KalmanFit::Helix::cosPhi0(), del4MXDelA(), KalmanFit::Helix::del4MXDelA(), delXDelA(), KalmanFit::Helix::delXDelA(), momentum(), KalmanFit::Helix::momentum(), operator=(), KalmanFit::Helix::operator=(), updateCache(), KalmanFit::Helix::updateCache(), x(), and KalmanFit::Helix::x().

HepSymMatrix Helix::m_Ea [private]

Definition at line 178 of file Helix.h.

Referenced by Ea(), KalmanFit::Helix::Ea(), ignoreErrorMatrix(), KalmanFit::Helix::ignoreErrorMatrix(), momentum(), KalmanFit::Helix::momentum(), operator=(), KalmanFit::Helix::operator=(), pivot(), KalmanFit::Helix::pivot(), set(), KalmanFit::Helix::set(), x(), and KalmanFit::Helix::x().

bool Helix::m_matrixValid [private]

Definition at line 179 of file Helix.h.

Referenced by ignoreErrorMatrix(), KalmanFit::Helix::ignoreErrorMatrix(), momentum(), KalmanFit::Helix::momentum(), operator=(), KalmanFit::Helix::operator=(), pivot(), KalmanFit::Helix::pivot(), set(), KalmanFit::Helix::set(), x(), and KalmanFit::Helix::x().

HepPoint3D Helix::m_pivot [private]

Definition at line 176 of file Helix.h.

Referenced by momentum(), KalmanFit::Helix::momentum(), operator=(), KalmanFit::Helix::operator=(), pivot(), KalmanFit::Helix::pivot(), set(), KalmanFit::Helix::set(), updateCache(), KalmanFit::Helix::updateCache(), x(), and KalmanFit::Helix::x().

IMagneticFieldSvc* Helix::m_pmgnIMF [protected]

Definition at line 162 of file Helix.h.

Referenced by Helix().

double Helix::m_pt [private]

Definition at line 185 of file Helix.h.

Referenced by momentum(), KalmanFit::Helix::momentum(), operator=(), KalmanFit::Helix::operator=(), pt(), updateCache(), and KalmanFit::Helix::updateCache().

double Helix::m_r [private]

Definition at line 186 of file Helix.h.

Referenced by curv(), KalmanFit::Helix::curv(), del4MXDelA(), KalmanFit::Helix::del4MXDelA(), delApDelA(), KalmanFit::Helix::delApDelA(), delXDelA(), KalmanFit::Helix::delXDelA(), momentum(), KalmanFit::Helix::momentum(), operator=(), KalmanFit::Helix::operator=(), pivot(), KalmanFit::Helix::pivot(), radius(), KalmanFit::Helix::radius(), updateCache(), KalmanFit::Helix::updateCache(), x(), and KalmanFit::Helix::x().

double Helix::m_sp [private]

Definition at line 184 of file Helix.h.

Referenced by del4MXDelA(), KalmanFit::Helix::del4MXDelA(), delXDelA(), KalmanFit::Helix::delXDelA(), momentum(), KalmanFit::Helix::momentum(), operator=(), KalmanFit::Helix::operator=(), sinPhi0(), KalmanFit::Helix::sinPhi0(), updateCache(), KalmanFit::Helix::updateCache(), x(), and KalmanFit::Helix::x().


Generated on Tue Nov 29 23:19:39 2016 for BOSS_7.0.2 by  doxygen 1.4.7