MdcxHel Class Reference

#include <MdcxHel.h>

Inheritance diagram for MdcxHel:

MdcxFittedHel MdcxSeg List of all members.

Public Member Functions

 MdcxHel ()
 MdcxHel (double D0, double Phi0, double Omega, double Z0=0.0, double Tanl=0.0, double T0=0.0, int Code=11111, int Mode=0, double X=0.0, double Y=0.0)
virtual ~MdcxHel ()
double D0 () const
double Phi0 () const
double Omega () const
double Z0 () const
double Tanl () const
double X () const
double Xref () const
double Y () const
double Yref () const
double T0 () const
double CosPhi0 () const
double SinPhi0 () const
double Doca_Len () const
double Doca_FLen () const
double Doca_Tof () const
double Doca_Zh () const
int Doca_Samb () const
int Doca_Wamb () const
double Doca_Eang () const
double Omin () const
int Mode () const
int Code () const
int Nfree () const
int Ominfl () const
int Qd0 () const
int Qphi0 () const
int Qomega () const
int Qz0 () const
int Qtanl () const
int Qt0 () const
double Xc () const
double Yc () const
double X0 () const
double Y0 () const
double Xh (double l) const
double Yh (double l) const
double Zh (double l) const
double Px (double l=0.0) const
double Py (double l=0.0) const
double Pz (double l=0.0) const
double Ptot (double l=0.0) const
double Lmax () const
double Doca (double WX, double WY, double WZ, double X, double Y, double Z=0.0)
double Doca (const MdcxHit &h)
void print () const
void flip ()
std::vector< float > derivatives (const MdcxHit &h)
void SetRef (double x, double y)
void SetMode (int n)
void SetD0 (int n)
void SetPhi0 (int n)
void SetOmega (int n)
void SetZ0 (int n)
void SetTanl (int n)
void SetT0 (int n)
void SetTurnFlag (const int &i)
int GetTurnFlag () const
MdcxHeloperator= (const MdcxHel &)

Protected Member Functions

double CalcPhi (double xf, double yf, double xl, double yl) const
int deltaq (int i, int j) const
void decode (const int i, int &i1, int &i2, int &i3, int &i4, int &i5, int &i6, int &n)
void copy (const MdcxHel &hel)

Protected Attributes

double d0
double phi0
double omega
double z0
double tanl
double t0
double xref
double yref
double cphi0
double sphi0
double x0
double y0
double xc
double yc
int code
int mode
int qd0
int qphi0
int qomega
int qz0
int qtanl
int qt0
int nfree
int ominfl
int turnflag
double omin
double len
double phi
double xh
double yh
double zh
double vx
double vy
double vz
double cosl
double sinl
double f0
double tx
double ty
double tz
Hep3Vector wvec
Hep3Vector tvec
Hep3Vector vhat
Hep3Vector dvec
int samb
int wamb
double eang

Detailed Description

Definition at line 37 of file MdcxHel.h.


Constructor & Destructor Documentation

MdcxHel::MdcxHel (  ) 

Definition at line 32 of file MdcxHel.cxx.

00032 { }

MdcxHel::MdcxHel ( double  D0,
double  Phi0,
double  Omega,
double  Z0 = 0.0,
double  Tanl = 0.0,
double  T0 = 0.0,
int  Code = 11111,
int  Mode = 0,
double  X = 0.0,
double  Y = 0.0 
)

Definition at line 35 of file MdcxHel.cxx.

References cos(), M_PI, phi0, sin(), and X0.

00036                                                     : 
00037 d0(D0), phi0(Phi0), omega(Omega), z0(Z0), tanl(Tanl), t0(T0),
00038 xref(X), yref(Y), code(Code), mode(Mode), omin(0.000005) { 
00039 
00040   ominfl = (fabs(omega) < omin) ? 0 : 1;
00041   double m_2pi = 2.0*M_PI;
00042   if (phi0 > M_PI)  phi0 -= m_2pi;
00043   if (phi0 < -M_PI) phi0 += m_2pi;
00044   cphi0 = cos(phi0);
00045   sphi0 = sin(phi0);
00046   x0 = X0();
00047   y0 = Y0();
00048   xc = Xc();
00049   yc = Yc();
00050   decode(code, qd0, qphi0, qomega, qz0, qtanl, qt0, nfree);
00051   turnflag = 1;  
00052   //std::cout << "MdcxHel::MdcxHel()  -> (x0, y0) (" << x0 << ", " << y0 << ")" << std::endl;
00053 }//endof MdcxHel

MdcxHel::~MdcxHel (  )  [virtual]

Definition at line 55 of file MdcxHel.cxx.

00055 { }


Member Function Documentation

double MdcxHel::CalcPhi ( double  xf,
double  yf,
double  xl,
double  yl 
) const [inline, protected]

Definition at line 167 of file MdcxHel.h.

References M_PI.

00167                                                                {
00168    double phit=atan2(yl-yf,xl-xf); return phit<0?phit+2*M_PI:phit;
00169   }//endof CalcPhi

int MdcxHel::Code (  )  const [inline]

Definition at line 74 of file MdcxHel.h.

References code.

Referenced by copy(), and MdcxFindTracks::TakeToOrigin().

00074 {return code;}

void MdcxHel::copy ( const MdcxHel hel  )  [protected]

Definition at line 219 of file MdcxHel.cxx.

References Code(), code, CosPhi0(), cphi0, D0(), d0, GetTurnFlag(), Mode(), mode, Nfree(), nfree, Omega(), omega, Omin(), omin, Ominfl(), ominfl, Phi0(), phi0, Qd0(), qd0, Qomega(), qomega, Qphi0(), qphi0, Qt0(), qt0, Qtanl(), qtanl, Qz0(), qz0, SinPhi0(), sphi0, T0(), t0, Tanl(), tanl, turnflag, X0(), x0, Xc(), xc, Xref(), xref, Y0(), y0, Yc(), yc, Yref(), yref, Z0(), and z0.

Referenced by MdcxFittedHel::Grow(), operator=(), and MdcxFittedHel::operator=().

00220 {
00221   //FIXME
00222   omega=rhs.Omega(); phi0=rhs.Phi0(); d0=rhs.D0(); t0=rhs.T0();
00223   tanl=rhs.Tanl(); z0=rhs.Z0();
00224   cphi0=rhs.CosPhi0(); sphi0=rhs.SinPhi0();
00225   x0=rhs.X0(); y0=rhs.Y0(); xc=rhs.Xc(); yc=rhs.Yc();
00226   xref=rhs.Xref(); yref=rhs.Yref();
00227   qomega=rhs.Qomega(); qphi0=rhs.Qphi0(); qd0=rhs.Qd0(); qt0=rhs.Qt0();
00228   qtanl=rhs.Qtanl(); qz0=rhs.Qz0();
00229   mode=rhs.Mode(); nfree=rhs.Nfree();
00230   code=rhs.Code(); ominfl=rhs.Ominfl(); omin=rhs.Omin();
00231   turnflag=rhs.GetTurnFlag();
00232 }//endof copy

double MdcxHel::CosPhi0 (  )  const [inline]

Definition at line 63 of file MdcxHel.h.

References cphi0.

Referenced by copy(), and MdcxFindTracks::TakeToOrigin().

00063 {return cphi0;}

double MdcxHel::D0 (  )  const [inline]

Definition at line 53 of file MdcxHel.h.

References d0.

Referenced by copy(), MdcxTrackFinder::FitMdcxTrack(), MdcxAddHits::GetAssociates(), MdcxFindTracks::process(), and MdcxFindTracks::TakeToOrigin().

00053 {return d0;}

void MdcxHel::decode ( const int  i,
int &  i1,
int &  i2,
int &  i3,
int &  i4,
int &  i5,
int &  i6,
int &  n 
) [protected]

Definition at line 199 of file MdcxHel.cxx.

References subSeperate::temp.

00201 {  
00202   int temp = code;
00203   temp=temp/1000000; temp=code-1000000*temp;
00204   i6=temp/100000;    temp=temp-100000*i6;
00205   i5=temp/10000;     temp=temp-10000*i5;
00206   i4=temp/1000;      temp=temp-1000*i4;
00207   i3=temp/100;       temp=temp-100*i3;
00208   i2=temp/10;        i1=temp-10*i2;
00209   n = 0;
00210   if(i6 == 1) n++; else i6 = 0;
00211   if(i5 == 1) n++; else i5 = 0;
00212   if(i4 == 1) n++; else i4 = 0;
00213   if(i3 == 1) n++; else i3 = 0;
00214   if(i2 == 1) n++; else i2 = 0;
00215   if(i1 == 1) n++; else i1 = 0;
00216 }//endof decode

int MdcxHel::deltaq ( int  i,
int  j 
) const [inline, protected]

Definition at line 170 of file MdcxHel.h.

Referenced by SetD0(), SetOmega(), SetPhi0(), SetT0(), SetTanl(), and SetZ0().

00170 { return i==j?0:i==0?1:-1; } //integer comparisons 

std::vector< float > MdcxHel::derivatives ( const MdcxHit h  ) 

Definition at line 324 of file MdcxHel.cxx.

References cos(), cosl, cphi0, d0, Doca(), f0, len, mode, nfree, omega, ominfl, phi, qd0, qomega, qphi0, qt0, qtanl, qz0, sin(), sphi0, subSeperate::temp, tx, ty, MdcxHit::v(), vx, vy, vz, x0, xh, y0, and yh.

Referenced by MdcxHit::derivatives(), and MdcxFittedHel::DoFit().

00325 {
00326   double doca = Doca(hit);
00327   std::vector<float> temp(nfree+1);
00328   temp[0] = doca;
00329   double fac = 1.0;
00330   if((mode==0) && (doca<0.0)) fac = -fac;
00331   if(mode == 0) temp[0] = fabs(temp[0]);
00332 
00333   int bump = 0;
00334   if (qd0)   temp[++bump] = (-vx*sphi0 + vy*cphi0) * fac;
00335   if (qphi0) {
00336     //double dddp0=-(yh-y0)*vx+(xh-x0)*vy;
00337     double dddp0 = -(yh-y0+f0*ty)*vx + (xh-x0+f0*tx)*vy;
00338     dddp0 *= (1.0 + d0*omega);
00339     temp[++bump] = dddp0*fac;
00340   }
00341   if (qomega) {
00342     double dddom;
00343     if (ominfl) {
00344       dddom = ((len*cos(phi)-xh+x0)*vx + (len*sin(phi)-yh+y0)*vy)/omega;
00345       dddom += f0*len*cosl*(-sin(phi)*vx+cos(phi)*vy);
00346     } else {
00347       dddom = 0.5*len*len*(-vx*sphi0+vy*cphi0);
00348     }
00349     temp[++bump] = dddom * fac;
00350   }
00351   if (qz0)   temp[++bump] = vz * fac;
00352   if (qtanl) temp[++bump] = (vz *len) * fac;
00353   if (qt0)   temp[++bump] = -hit.v();
00354   return temp;
00355 }//endof derivatives

double MdcxHel::Doca ( const MdcxHit h  ) 

Definition at line 234 of file MdcxHel.cxx.

References Doca(), MdcxHit::wx(), MdcxHit::wy(), MdcxHit::wz(), MdcxHit::x(), and MdcxHit::y().

00234                                      { 
00235   //std::cout<< __FILE__ << "   " << __LINE__ << "  hit("<<h.Layer()<<","<<h.WireNo()<<")";
00236   return Doca( h.wx(), h.wy(), h.wz(), h.x(), h.y() ); 
00237 }

double MdcxHel::Doca ( double  WX,
double  WY,
double  WZ,
double  X,
double  Y,
double  Z = 0.0 
)

Definition at line 239 of file MdcxHel.cxx.

References cos(), cosl, cphi0, dvec, eang, f0, f1, f2, len, M_PI, MdcxParameters::maxMdcZLen, omega, ominfl, phi, phi0, Constants::pi, samb, sin(), sinl, sphi0, tanl, turnflag, tvec, tx, ty, tz, vhat, vx, vy, vz, wamb, wvec, xc, Xh(), xh, xref, yc, Yh(), yh, yref, z0, Zh(), and zh.

Referenced by MdcxHit::d(), derivatives(), Doca(), MdcxFindTracks::drophits(), MdcxFindTracks::process(), MdcxHit::pull(), MdcxHit::residual(), and MdcxFindTracks::resout().

00241 {
00242   double m_2pi = 2.0*M_PI;
00243   // describe wire
00244   //cout << " In Doca, xi = " << xi << " yi = " << yi << " zi = " << zi <<endl;
00245   Hep3Vector ivec(xi, yi, zi); 
00246   wvec = Hep3Vector(wx, wy, wz);
00247   //cout << " In Doca, wx = " << wx << " wy = " << wy << " wz = " << wz <<endl;
00248   //  calculate len to doca
00249   double zd, xd = xi, yd = yi;
00250   // cout << " In Doca, start xd = " << xd << " yd = " << yd << endl;
00251   double lnew,t1,t2,dphi,dlen=1000.0;
00252   len = 0.0;
00253   int itry = 2;
00254   // int segflg=0; if ((code==111)&&(z0==0.0)&&(tanl==0.0))segflg=1;
00255   // int superseg=0; if ((code==11111)&&(xref!=0.0)&&(yref!=0.0))superseg=1;
00256   double circut, circum = 10000.; 
00257   if (ominfl) circum = m_2pi/fabs(omega);
00258   circut = 0.50 * circum;
00259   while (itry) {
00260     if (ominfl) {
00262       t1 = -xc + xd; t2 = yc - yd; phi = atan2(t1, t2); 
00263       if (omega < 0.0) phi += M_PI;
00264       if (phi > M_PI)  phi -= m_2pi;
00265       dphi = phi - phi0;
00266       if (omega < 0.0){
00267         if (dphi > 0.0)    dphi -= m_2pi;
00268         if (dphi < -m_2pi) dphi += m_2pi;
00269       }else{
00270         if (dphi < 0.0)   dphi += m_2pi;
00271         if (dphi > m_2pi) dphi -= m_2pi;
00272       }
00273       lnew = dphi/omega; 
00274       //   if ((lnew>circut)&&(segflg))lnew-=circum; 
00275       //   if ((lnew>circut)&&(superseg))lnew-=circum; 
00276       if ((lnew>circut)&&(turnflag)) lnew -= circum;  //FIXME attention
00277       
00278       zh = Zh(lnew); 
00279       xd=xi+(zh-zi)*wx/wz; yd=yi+(zh-zi)*wy/wz; zd=zh;
00280       //   cout << " In Doca, xd = " << xd << " yd = " << yd << " zh = " << zh;
00281       //   cout << " lnew = " << lnew << endl;
00282       dlen=fabs(lnew-len); len=lnew; 
00283       //   if (segflg)break; 
00284       //std::cout<< __FILE__ << __LINE__<<" Doca()  dlen " << dlen<< " zh  "<<zh<<" >?"
00285         //<<MdcxParameters::maxMdcZLen<<std::endl;
00286       if (fabs(zh) > MdcxParameters::maxMdcZLen) break;  //FIXME attention
00287       if ( (0.0==wx) && (0.0==wy) )break; if (dlen < 0.000001)break; itry--;
00288     } else {
00289       len = (xi-xref)*cphi0 + (yi-yref)*sphi0;
00290       zh  = z0 + tanl*len;
00291       phi = phi0;
00292       break;
00293     }
00294   }
00295   //  Hep3Vector Dvec(xd,yd,zd);
00296   xh = Xh(len); yh = Yh(len); 
00297   Hep3Vector hvec(xh, yh, zh);
00298 //   cout << " In Doca, xh = " << xh << " yh = " << yh << " zh = " << zh << " len=" << len << " om " << omega << endl;
00299   double lamb = atan(tanl);
00300   cosl = cos(lamb); sinl = sin(lamb);
00301   tx = cosl*cos(phi); ty = cosl*sin(phi); tz = sinl; 
00302   tvec = Hep3Vector(tx, ty, tz); 
00303   Hep3Vector vvec = wvec.cross(tvec); 
00304   vhat = vvec.unit(); vx = vhat.x(); vy = vhat.y(); vz = vhat.z();
00305   // cout << " In Doca, vx = " << vx << " vy = " << vy << " vz = " << vz << endl;
00306   dvec = ivec - hvec;
00307   double doca = dvec*vhat;
00308   // cout << " doca = " << doca << endl;
00309   double f1 = dvec*tvec; double f2 = wvec*tvec; double f3 = dvec*wvec;
00310   f0 = (f1 - f2*f3) / (1.0 - f2*f2); 
00311   samb = (doca > 0.0) ? -1 : +1;
00312   double wirephi = atan2(yd, xd);
00313   eang = BesAngle(phi-wirephi);
00314   wamb = (fabs(eang) < Constants::pi/2) ? samb : -samb;
00315   //std::cout<< __FILE__ << __LINE__<<" Doca()  dlen " << dlen<< " zh  "<<zh<<" >?"
00316         //<<MdcxParameters::maxMdcZLen<<std::endl;
00317   if (fabs(zh) > MdcxParameters::maxMdcZLen) doca = 1000.0; 
00318   //if(doca == 1000.0)    cout << " In Doca, zh = " << zh << " len=" << len << " om " << omega <<" "<< ominfl<< 
00319    //" z0 " << z0 << "tanl " << tanl <<endl;
00320   //cout << "  doca = " << doca << endl;
00321   return doca;
00322 }//endof Doca

double MdcxHel::Doca_Eang (  )  const [inline]

Definition at line 71 of file MdcxHel.h.

References eang.

Referenced by MdcxHit::d(), MdcxHit::derivatives(), MdcxHit::pull(), and MdcxHit::residual().

00071 {return eang;}

double MdcxHel::Doca_FLen (  )  const [inline]

Definition at line 66 of file MdcxHel.h.

References len, and tanl.

00066 {return len*sqrt(1.0+tanl*tanl);}

double MdcxHel::Doca_Len (  )  const [inline]

Definition at line 65 of file MdcxHel.h.

References len.

Referenced by MdcxAddHits::GetAssociates(), and MdcxFindTracks::process().

00065 {return len;}

int MdcxHel::Doca_Samb (  )  const [inline]

Definition at line 69 of file MdcxHel.h.

References samb.

00069 {return samb;}

double MdcxHel::Doca_Tof (  )  const [inline]

Definition at line 67 of file MdcxHel.h.

References MdcxParameters::c, len, and tanl.

Referenced by MdcxHit::d(), MdcxHit::derivatives(), MdcxHit::pull(), MdcxHit::residual(), and MdcxFindTracks::resout().

00067 {return len*sqrt(1.0+tanl*tanl)/MdcxParameters::c;}

int MdcxHel::Doca_Wamb (  )  const [inline]

Definition at line 70 of file MdcxHel.h.

References wamb.

Referenced by MdcxHit::d(), MdcxHit::derivatives(), MdcxHit::pull(), and MdcxHit::residual().

00070 {return wamb;}

double MdcxHel::Doca_Zh (  )  const [inline]

Definition at line 68 of file MdcxHel.h.

References zh.

Referenced by MdcxHit::d(), MdcxHit::derivatives(), MdcxHit::pull(), MdcxHit::residual(), and MdcxFindTracks::resout().

00068 {return zh;}

void MdcxHel::flip (  ) 

Definition at line 375 of file MdcxHel.cxx.

References cos(), cphi0, d0, M_PI, omega, ominfl, phi0, sin(), sphi0, tanl, X0(), x0, Y0(), y0, z0, and Zh().

Referenced by MdcxFindTracks::process().

00375                    {
00376   double m_2pi = 2.0*M_PI;
00377   if (ominfl) {
00378     if ( (fabs(d0) + 2.0/fabs(omega)) > 80.0 ) return;
00379     double lturn = m_2pi/fabs(omega);
00380     double zturn = Zh(lturn);
00381     //  cout << "z0 " << z0 << " zturn " << zturn << endl;
00382     if (fabs(zturn) < fabs(z0)) {
00383       z0 = zturn;
00384       tanl = -tanl;
00385       omega = -omega;
00386       d0 = -d0; 
00387       phi0 = phi0 - M_PI;
00388       if (phi0 < -M_PI) phi0 += m_2pi;
00389       cphi0 = cos(phi0);
00390       sphi0 = sin(phi0);
00391       x0 = X0();
00392       y0 = Y0();
00393     } 
00394   }
00395 }//endof flip

int MdcxHel::GetTurnFlag (  )  const [inline]

Definition at line 116 of file MdcxHel.h.

References turnflag.

Referenced by copy().

00116 {return turnflag;}

double MdcxHel::Lmax (  )  const

Definition at line 142 of file MdcxHel.cxx.

References d0, dmax, M_PI, MdcxParameters::maxMdcRadius, MdcxParameters::maxTrkLength, omega, and ominfl.

Referenced by MdcxAddHits::GetAssociates(), and MdcxFindTracks::process().

00142                            {
00143   double lmax = MdcxParameters::maxTrkLength;
00144   if (ominfl) {
00145     double rmax = 1.0/fabs(omega);
00146     double dmax = fabs(d0) + 2.0*rmax;
00147     if (dmax > MdcxParameters::maxMdcRadius) lmax = M_PI*rmax;  
00148   }
00149   return lmax;
00150 }//endof Lmax

int MdcxHel::Mode (  )  const [inline]

Definition at line 73 of file MdcxHel.h.

References mode.

Referenced by copy(), MdcxHit::pull(), and MdcxHit::residual().

00073 {return mode;}

int MdcxHel::Nfree (  )  const [inline]

Definition at line 75 of file MdcxHel.h.

References nfree.

Referenced by copy().

00075 {return nfree;}

double MdcxHel::Omega (  )  const [inline]

Definition at line 55 of file MdcxHel.h.

References omega.

Referenced by copy(), MdcxTrackFinder::FitMdcxTrack(), MdcxAddHits::GetAssociates(), MdcxFindTracks::process(), and MdcxFindTracks::TakeToOrigin().

00055 {return omega;}

double MdcxHel::Omin (  )  const [inline]

Definition at line 72 of file MdcxHel.h.

References omin.

Referenced by copy().

00072 {return omin;}

int MdcxHel::Ominfl (  )  const [inline]

Definition at line 76 of file MdcxHel.h.

References ominfl.

Referenced by copy(), and MdcxAddHits::GetAssociates().

00076 {return ominfl;}

MdcxHel & MdcxHel::operator= ( const MdcxHel  ) 

Reimplemented in MdcxFittedHel.

Definition at line 193 of file MdcxHel.cxx.

References copy().

00193                                               {
00194   copy(rhs);
00195   return *this;
00196 }

double MdcxHel::Phi0 (  )  const [inline]

Definition at line 54 of file MdcxHel.h.

References phi0.

Referenced by copy(), MdcxTrackFinder::FitMdcxTrack(), MdcxAddHits::GetAssociates(), MdcxFindTracks::process(), and MdcxFindTracks::TakeToOrigin().

00054 {return phi0;}

void MdcxHel::print (  )  const

Definition at line 357 of file MdcxHel.cxx.

References code, d0, mode, nfree, omega, ominfl, phi0, t0, tanl, x0, xc, xref, y0, yc, yref, and z0.

Referenced by MdcxFittedHel::DoFit(), MdcxFittedHel::FitPrint(), MdcxAddHits::GetAssociates(), MdcxFittedHel::IterateFit(), and MdcxFindTracks::process().

00357                           {
00358   cout << "MdcxHel(";
00359   cout << d0<<",";
00360   cout << phi0<<",";
00361   cout << omega<<",";
00362   cout << z0<<",";
00363   cout << tanl<<")"<<endl; 
00364   cout << " t0 = " << t0 ; 
00365   cout << " nfree = " << nfree ;
00366   cout << " (x0,y0) " << x0<<","<<y0;
00367   cout << " (xc,yc) " << xc<<","<<yc;
00368   cout << " (xref,yref) " << xref<<","<<yref;
00369   cout << " code = " << code; 
00370   cout << " mode = " << mode;
00371   cout << " ominfl = " << ominfl;
00372   cout << "  " << endl;
00373 }//endof print

double MdcxHel::Ptot ( double  l = 0.0  )  const

Definition at line 134 of file MdcxHel.cxx.

References omega, ominfl, and tanl.

00134                                  {
00135   if(ominfl) {
00136     return 0.003*sqrt(1.0+tanl*tanl)/fabs(omega);
00137   } else {
00138     return 1000.0*sqrt(1.0+tanl*tanl);
00139   }//ominfl
00140 }//endof Ptot

double MdcxHel::Px ( double  l = 0.0  )  const

Definition at line 107 of file MdcxHel.cxx.

References cos(), cphi0, omega, ominfl, and phi0.

00107                                  {
00108   if(ominfl) {
00109     double phit = phi0 + omega*l;
00110     return 0.003*cos(phit)/fabs(omega);  
00111   } else {
00112     return 1000.0*cphi0;
00113   }//ominfl
00114 }//endof Px

double MdcxHel::Py ( double  l = 0.0  )  const

Definition at line 116 of file MdcxHel.cxx.

References omega, ominfl, phi0, sin(), and sphi0.

00116                                  {
00117   if(ominfl) {
00118     double phit = phi0+omega*l;
00119     return 0.003*sin(phit)/fabs(omega);
00120   } else {
00121     return 1000.0*sphi0;
00122   }//ominfl
00123 }//endof Py

double MdcxHel::Pz ( double  l = 0.0  )  const

Definition at line 125 of file MdcxHel.cxx.

References omega, ominfl, and tanl.

00125                                  {
00126   if(ominfl) {
00127     return 0.003*tanl/fabs(omega);
00128   }
00129   else{
00130     return 1000.0*tanl;
00131   }//ominfl
00132 }//endof Pz

int MdcxHel::Qd0 (  )  const [inline]

Definition at line 77 of file MdcxHel.h.

References qd0.

Referenced by copy().

00077 {return qd0;}

int MdcxHel::Qomega (  )  const [inline]

Definition at line 79 of file MdcxHel.h.

References qomega.

Referenced by copy().

00079 {return qomega;}

int MdcxHel::Qphi0 (  )  const [inline]

Definition at line 78 of file MdcxHel.h.

References qphi0.

Referenced by copy().

00078 {return qphi0;}

int MdcxHel::Qt0 (  )  const [inline]

Definition at line 82 of file MdcxHel.h.

References qt0.

Referenced by copy().

00082 {return qt0;}

int MdcxHel::Qtanl (  )  const [inline]

Definition at line 81 of file MdcxHel.h.

References qtanl.

Referenced by copy().

00081 {return qtanl;}

int MdcxHel::Qz0 (  )  const [inline]

Definition at line 80 of file MdcxHel.h.

References qz0.

Referenced by copy().

00080 {return qz0;}

void MdcxHel::SetD0 ( int  n  ) 

Definition at line 171 of file MdcxHel.cxx.

References code, deltaq(), nfree, and qd0.

00171                            {
00172   nfree = nfree + deltaq(qd0, Qd0);
00173   code  = code + deltaq(qd0, Qd0);
00174   qd0 = Qd0;
00175 }

void MdcxHel::SetMode ( int  n  ) 

Definition at line 154 of file MdcxHel.cxx.

References mode.

00154 { mode = n; }

void MdcxHel::SetOmega ( int  n  ) 

Definition at line 161 of file MdcxHel.cxx.

References code, deltaq(), nfree, and qomega.

00161                                  {
00162   nfree = nfree + deltaq(qomega, Qomega);
00163   code  = code + deltaq(qomega, Qomega)*100;
00164   qomega = Qomega;
00165 }

void MdcxHel::SetPhi0 ( int  n  ) 

Definition at line 166 of file MdcxHel.cxx.

References code, deltaq(), nfree, and qphi0.

00166                                {
00167   nfree = nfree + deltaq(qphi0, Qphi0);
00168   code  = code + deltaq(qphi0, Qphi0)*10;
00169   qphi0 = Qphi0;
00170 }

void MdcxHel::SetRef ( double  x,
double  y 
)

Definition at line 156 of file MdcxHel.cxx.

References xref, and yref.

00156                                        {
00157   xref = x;
00158   yref = y;
00159 }

void MdcxHel::SetT0 ( int  n  ) 

Definition at line 186 of file MdcxHel.cxx.

References code, deltaq(), nfree, and qt0.

00186                            {
00187   nfree = nfree + deltaq(qt0, Qt0);
00188   code = code + deltaq(qt0, Qt0)*100000;
00189   qt0 = Qt0;
00190 }

void MdcxHel::SetTanl ( int  n  ) 

Definition at line 176 of file MdcxHel.cxx.

References code, deltaq(), nfree, and qtanl.

00176                                {
00177   nfree = nfree + deltaq(qtanl, Qtanl);
00178   code  = code + deltaq(qtanl, Qtanl)*10000;
00179   qtanl = Qtanl;
00180 }

void MdcxHel::SetTurnFlag ( const int &  i  )  [inline]

Definition at line 115 of file MdcxHel.h.

References turnflag.

Referenced by MdcxMergeDups::MdcxMergeDups(), and MdcxFindSegs::trial().

00115 {turnflag=i;}

void MdcxHel::SetZ0 ( int  n  ) 

Definition at line 181 of file MdcxHel.cxx.

References code, deltaq(), nfree, and qz0.

00181                            {
00182   nfree = nfree + deltaq(qz0, Qz0);
00183   code = code + deltaq(qz0, Qz0)*1000;
00184   qz0 = Qz0;
00185 }

double MdcxHel::SinPhi0 (  )  const [inline]

Definition at line 64 of file MdcxHel.h.

References sphi0.

Referenced by copy(), and MdcxFindTracks::TakeToOrigin().

00064 {return sphi0;}

double MdcxHel::T0 (  )  const [inline]

Definition at line 62 of file MdcxHel.h.

References t0.

Referenced by copy(), MdcxHit::d(), MdcxHit::derivatives(), MdcxHit::pull(), MdcxHit::residual(), and MdcxFindTracks::TakeToOrigin().

00062 {return t0;}

double MdcxHel::Tanl (  )  const [inline]

Definition at line 57 of file MdcxHel.h.

References tanl.

Referenced by copy(), MdcxTrackFinder::FitMdcxTrack(), and MdcxFindTracks::TakeToOrigin().

00057 {return tanl;}

double MdcxHel::X (  )  const [inline]

Definition at line 58 of file MdcxHel.h.

References xref.

00058 {return xref;}

double MdcxHel::X0 (  )  const

Definition at line 77 of file MdcxHel.cxx.

References d0, sphi0, and xref.

Referenced by copy(), MdcxFittedHel::DoFit(), flip(), and MdcxFindTracks::process().

00077                        {
00078   return (xref - sphi0*d0);
00079 }//endof X0

double MdcxHel::Xc (  )  const

Definition at line 59 of file MdcxHel.cxx.

References omega, ominfl, sphi0, and x0.

Referenced by copy(), MdcxFittedHel::DoFit(), MdcxAddHits::GetAssociates(), and MdcxFindTracks::process().

00059                         {
00060   if(ominfl) {
00061     //return (X0() - sphi0/omega);
00062     return (x0 - sphi0/omega);
00063   } else {
00064     return 999999999.9;
00065   }//(ominfl)
00066 }//endof Xc

double MdcxHel::Xh ( double  l  )  const

Definition at line 85 of file MdcxHel.cxx.

References cphi0, omega, ominfl, phi0, sin(), sphi0, x0, and xc.

Referenced by Doca(), and MdcxAddHits::GetAssociates().

00085                                {
00086   if(ominfl){
00087     double phit=phi0+omega*l;
00088     return (xc+sin(phit)/omega);
00089   }else{
00090     return (x0+cphi0*l-0.5*l*l*omega*sphi0);
00091   }//ominfl
00092 }//endof Xh

double MdcxHel::Xref (  )  const [inline]

Definition at line 59 of file MdcxHel.h.

References xref.

Referenced by copy(), MdcxSeg::printSegAll(), and MdcxFindTracks::TakeToOrigin().

00059 {return xref;}

double MdcxHel::Y (  )  const [inline]

Definition at line 60 of file MdcxHel.h.

References yref.

00060 {return yref;}

double MdcxHel::Y0 (  )  const

Definition at line 81 of file MdcxHel.cxx.

References cphi0, d0, and yref.

Referenced by copy(), MdcxFittedHel::DoFit(), flip(), and MdcxFindTracks::process().

00081                        {
00082   return (yref + cphi0*d0);
00083 }//endof Y0

double MdcxHel::Yc (  )  const

Definition at line 68 of file MdcxHel.cxx.

References cphi0, omega, ominfl, and y0.

Referenced by copy(), MdcxFittedHel::DoFit(), MdcxAddHits::GetAssociates(), and MdcxFindTracks::process().

00068                         {
00069   if(ominfl) {
00070     //return (Y0()+cphi0/omega);
00071     return (y0 + cphi0/omega);
00072   } else {
00073     return 999999999.9;
00074   }//(ominfl)
00075 }//endof Yc

double MdcxHel::Yh ( double  l  )  const

Definition at line 94 of file MdcxHel.cxx.

References cos(), cphi0, omega, ominfl, phi0, sphi0, y0, and yc.

Referenced by Doca(), and MdcxAddHits::GetAssociates().

00094                                {
00095   if(ominfl){
00096     double phit=phi0+omega*l;
00097     return (yc-cos(phit)/omega);
00098   }else{
00099     return (y0+sphi0*l+0.5*l*l*omega*cphi0);
00100   }//ominfl
00101 }//endof Yh

double MdcxHel::Yref (  )  const [inline]

Definition at line 61 of file MdcxHel.h.

References yref.

Referenced by copy(), MdcxSeg::printSegAll(), and MdcxFindTracks::TakeToOrigin().

00061 {return yref;}

double MdcxHel::Z0 (  )  const [inline]

Definition at line 56 of file MdcxHel.h.

References z0.

Referenced by copy(), MdcxTrackFinder::FitMdcxTrack(), and MdcxFindTracks::TakeToOrigin().

00056 {return z0;}

double MdcxHel::Zh ( double  l  )  const

Definition at line 103 of file MdcxHel.cxx.

References tanl, and z0.

Referenced by Doca(), and flip().

00103                                {
00104   return (z0+tanl*l);
00105 }//endof Zh


Member Data Documentation

int MdcxHel::code [protected]

Definition at line 137 of file MdcxHel.h.

Referenced by Code(), copy(), print(), SetD0(), SetOmega(), SetPhi0(), SetT0(), SetTanl(), and SetZ0().

double MdcxHel::cosl [protected]

Definition at line 162 of file MdcxHel.h.

Referenced by derivatives(), and Doca().

double MdcxHel::cphi0 [protected]

Definition at line 131 of file MdcxHel.h.

Referenced by copy(), CosPhi0(), derivatives(), Doca(), MdcxFittedHel::DoFit(), flip(), Px(), Xh(), Y0(), Yc(), Yh(), and MdcxSeg::Yline_bbrrf().

double MdcxHel::d0 [protected]

Definition at line 123 of file MdcxHel.h.

Referenced by copy(), D0(), MdcxSeg::D0_sl_approx(), derivatives(), MdcxFittedHel::DoFit(), MdcxFittedHel::FitPrint(), flip(), Lmax(), print(), MdcxSeg::printSegAll(), X0(), MdcxSeg::Xline_bbrrf(), Y0(), and MdcxSeg::Yline_bbrrf().

Hep3Vector MdcxHel::dvec [protected]

Definition at line 163 of file MdcxHel.h.

Referenced by Doca().

double MdcxHel::eang [protected]

Definition at line 164 of file MdcxHel.h.

Referenced by Doca(), and Doca_Eang().

double MdcxHel::f0 [protected]

Definition at line 162 of file MdcxHel.h.

Referenced by derivatives(), and Doca().

double MdcxHel::len [protected]

Definition at line 162 of file MdcxHel.h.

Referenced by derivatives(), Doca(), Doca_FLen(), Doca_Len(), and Doca_Tof().

int MdcxHel::mode [protected]

Definition at line 145 of file MdcxHel.h.

Referenced by copy(), derivatives(), Mode(), print(), and SetMode().

int MdcxHel::nfree [protected]

Definition at line 153 of file MdcxHel.h.

Referenced by copy(), derivatives(), MdcxFittedHel::DoFit(), MdcxFittedHel::Grow(), MdcxFittedHel::IterateFit(), Nfree(), print(), SetD0(), SetOmega(), SetPhi0(), SetT0(), SetTanl(), and SetZ0().

double MdcxHel::omega [protected]

Definition at line 125 of file MdcxHel.h.

Referenced by copy(), MdcxSeg::D0_sl_approx(), derivatives(), Doca(), MdcxFittedHel::DoFit(), MdcxFittedHel::FitPrint(), flip(), Lmax(), Omega(), MdcxSeg::Phi0_sl_approx(), print(), MdcxSeg::printSegAll(), Ptot(), Px(), Py(), Pz(), Xc(), Xh(), Yc(), and Yh().

double MdcxHel::omin [protected]

Definition at line 158 of file MdcxHel.h.

Referenced by copy(), MdcxFittedHel::DoFit(), and Omin().

int MdcxHel::ominfl [protected]

Definition at line 154 of file MdcxHel.h.

Referenced by copy(), derivatives(), Doca(), MdcxFittedHel::DoFit(), flip(), Lmax(), Ominfl(), print(), Ptot(), Px(), Py(), Pz(), Xc(), Xh(), Yc(), and Yh().

double MdcxHel::phi [protected]

Definition at line 162 of file MdcxHel.h.

Referenced by derivatives(), and Doca().

double MdcxHel::phi0 [protected]

Definition at line 124 of file MdcxHel.h.

Referenced by copy(), Doca(), MdcxFittedHel::DoFit(), MdcxFittedHel::FitPrint(), flip(), Phi0(), MdcxSeg::Phi0_sl_approx(), print(), MdcxSeg::printSegAll(), Px(), Py(), Xh(), and Yh().

int MdcxHel::qd0 [protected]

Definition at line 147 of file MdcxHel.h.

Referenced by copy(), derivatives(), MdcxFittedHel::DoFit(), Qd0(), and SetD0().

int MdcxHel::qomega [protected]

Definition at line 149 of file MdcxHel.h.

Referenced by copy(), derivatives(), MdcxFittedHel::DoFit(), Qomega(), and SetOmega().

int MdcxHel::qphi0 [protected]

Definition at line 148 of file MdcxHel.h.

Referenced by copy(), derivatives(), MdcxFittedHel::DoFit(), Qphi0(), and SetPhi0().

int MdcxHel::qt0 [protected]

Definition at line 152 of file MdcxHel.h.

Referenced by copy(), derivatives(), MdcxFittedHel::DoFit(), Qt0(), and SetT0().

int MdcxHel::qtanl [protected]

Definition at line 151 of file MdcxHel.h.

Referenced by copy(), derivatives(), MdcxFittedHel::DoFit(), Qtanl(), and SetTanl().

int MdcxHel::qz0 [protected]

Definition at line 150 of file MdcxHel.h.

Referenced by copy(), derivatives(), MdcxFittedHel::DoFit(), Qz0(), and SetZ0().

int MdcxHel::samb [protected]

Definition at line 164 of file MdcxHel.h.

Referenced by Doca(), and Doca_Samb().

double MdcxHel::sinl [protected]

Definition at line 162 of file MdcxHel.h.

Referenced by Doca().

double MdcxHel::sphi0 [protected]

Definition at line 132 of file MdcxHel.h.

Referenced by copy(), derivatives(), Doca(), MdcxFittedHel::DoFit(), flip(), Py(), SinPhi0(), X0(), Xc(), Xh(), MdcxSeg::Xline_bbrrf(), and Yh().

double MdcxHel::t0 [protected]

Definition at line 128 of file MdcxHel.h.

Referenced by copy(), MdcxFittedHel::DoFit(), print(), and T0().

double MdcxHel::tanl [protected]

Definition at line 127 of file MdcxHel.h.

Referenced by copy(), Doca(), Doca_FLen(), Doca_Tof(), MdcxFittedHel::DoFit(), MdcxFittedHel::FitPrint(), flip(), print(), Ptot(), Pz(), Tanl(), and Zh().

int MdcxHel::turnflag [protected]

Definition at line 155 of file MdcxHel.h.

Referenced by copy(), Doca(), GetTurnFlag(), and SetTurnFlag().

Hep3Vector MdcxHel::tvec [protected]

Definition at line 163 of file MdcxHel.h.

Referenced by Doca().

double MdcxHel::tx [protected]

Definition at line 162 of file MdcxHel.h.

Referenced by derivatives(), and Doca().

double MdcxHel::ty [protected]

Definition at line 162 of file MdcxHel.h.

Referenced by derivatives(), and Doca().

double MdcxHel::tz [protected]

Definition at line 162 of file MdcxHel.h.

Referenced by Doca().

Hep3Vector MdcxHel::vhat [protected]

Definition at line 163 of file MdcxHel.h.

Referenced by Doca().

double MdcxHel::vx [protected]

Definition at line 162 of file MdcxHel.h.

Referenced by derivatives(), and Doca().

double MdcxHel::vy [protected]

Definition at line 162 of file MdcxHel.h.

Referenced by derivatives(), and Doca().

double MdcxHel::vz [protected]

Definition at line 162 of file MdcxHel.h.

Referenced by derivatives(), and Doca().

int MdcxHel::wamb [protected]

Definition at line 164 of file MdcxHel.h.

Referenced by Doca(), and Doca_Wamb().

Hep3Vector MdcxHel::wvec [protected]

Definition at line 163 of file MdcxHel.h.

Referenced by Doca().

double MdcxHel::x0 [protected]

Definition at line 133 of file MdcxHel.h.

Referenced by copy(), derivatives(), MdcxFittedHel::DoFit(), flip(), print(), Xc(), and Xh().

double MdcxHel::xc [protected]

Definition at line 135 of file MdcxHel.h.

Referenced by copy(), Doca(), MdcxFittedHel::DoFit(), print(), and Xh().

double MdcxHel::xh [protected]

Definition at line 162 of file MdcxHel.h.

Referenced by derivatives(), and Doca().

double MdcxHel::xref [protected]

Definition at line 129 of file MdcxHel.h.

Referenced by copy(), Doca(), print(), SetRef(), X(), X0(), MdcxSeg::Xline_bbrrf(), and Xref().

double MdcxHel::y0 [protected]

Definition at line 134 of file MdcxHel.h.

Referenced by copy(), derivatives(), MdcxFittedHel::DoFit(), flip(), print(), Yc(), and Yh().

double MdcxHel::yc [protected]

Definition at line 136 of file MdcxHel.h.

Referenced by copy(), Doca(), MdcxFittedHel::DoFit(), print(), and Yh().

double MdcxHel::yh [protected]

Definition at line 162 of file MdcxHel.h.

Referenced by derivatives(), and Doca().

double MdcxHel::yref [protected]

Definition at line 130 of file MdcxHel.h.

Referenced by copy(), Doca(), print(), SetRef(), Y(), Y0(), MdcxSeg::Yline_bbrrf(), and Yref().

double MdcxHel::z0 [protected]

Definition at line 126 of file MdcxHel.h.

Referenced by copy(), Doca(), MdcxFittedHel::DoFit(), MdcxFittedHel::FitPrint(), flip(), print(), Z0(), and Zh().

double MdcxHel::zh [protected]

Definition at line 162 of file MdcxHel.h.

Referenced by Doca(), and Doca_Zh().


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