00001 #include <stdlib.h>
00002 #include <math.h>
00003 #include "TError.h"
00004 #include "TArrayD.h"
00005 #include "TCanvas.h"
00006 #include "TGraph.h"
00007 #include "TSystem.h"
00008 #include "TMath.h"
00009 #if ROOT_VERSION_CODE < 331013
00010 #include "TCL.h"
00011 #else
00012 #include "TCernLib.h"
00013 #endif
00014 #include "TRandom.h"
00015 #include "TRandom2.h"
00016 #include "TError.h"
00017 #include "THelixTrack.h"
00018 #include "StMatrixD.hh"
00019 #include "TComplex.h"
00020 #include "TH1.h"
00021 #include <cassert>
00022
00023
00024 const TComplex Im(0,1);
00025
00026 inline static double dot(const TComplex &a,const TComplex &b)
00027 {return a.Re()*b.Re()+a.Im()*b.Im();}
00028
00029 inline static TComplex expOne(TComplex x)
00030 {
00031 double a = TComplex::Abs(x);
00032 if (a<0.01) {
00033 return 1.+x*((1/2.) + x*((1/6.)+ x*(1/24.)));
00034 } else {
00035 return (TComplex::Exp(x)-1.)/x;
00036 }
00037 }
00038
00039 inline static TComplex expOneD(TComplex x)
00040 {
00041 double a = TComplex::Abs(x);
00042 if (a<0.01) {
00043 return (1/2. + x*((1/6.)+ x*((1/24.)+x*(1/120.))));
00044 } else {
00045 return (TComplex::Exp(x)-1.-x)/(x*x);
00046 }
00047 }
00048
00049
00050 void TCEmx_t::Set(const double *err)
00051 { if (err) {memcpy(this,err,sizeof(*this));} else {Clear();}}
00052
00053
00054 void TCEmx_t::Move(const double F[3][3])
00055 {
00056 assert(mHH);
00057 double oErr[6];
00058 memcpy(oErr,Arr(),sizeof(oErr));
00059 TCL::trasat(F[0],oErr,Arr(),3,3);
00060 }
00061
00062 void TCEmx_t::Backward()
00063 {
00064 mHA*=-1; mAC*=-1;
00065 }
00066
00067 double TCEmx_t::Sign() const
00068 {
00069 const double *E = &mHH;
00070 double dia[3];
00071 for (int i=0,li=0;i<3;li+=++i) {
00072 dia[i] = E[li+i];
00073 if (dia[i]<0) return -(i+1)*11;
00074 for (int j=0;j<=i;j++) {
00075 double dis = dia[i]*dia[j]-E[li+j]*E[li+j];
00076 if (dis<0) return -(i+1+10*(j+1));
00077 } }
00078 return 0;
00079
00080 }
00081
00082
00083 void THEmx_t::Set(const double *errxy,const double *errz)
00084 {
00085 Clear();
00086 memcpy(&mHH,errxy,sizeof(mHH)*6);
00087 mZZ = errz[0]; mZL =errz[1]; mLL =errz[2];
00088 }
00089
00090 void THEmx_t::Set(const double *err)
00091 { if (err) {memcpy(this,err,sizeof(*this));} else {Clear();}}
00092
00093 void THEmx_t::Backward()
00094 {
00095 mHA*=-1; mAC*=-1; mHZ*=-1; mCZ*=-1; mAL*=-1; mZL*=-1;
00096 }
00097
00098 void THEmx_t::Move(const double F[5][5])
00099 {
00100 assert(mHH);
00101 double oErr[15];
00102 memcpy(oErr,Arr(),sizeof(oErr));
00103 TCL::trasat(F[0],oErr,Arr(),5,5);
00104 }
00105
00106 void THEmx_t::Print(const char *tit) const
00107 {
00108 static const char *N="HACZL";
00109 if (!tit) tit = "";
00110 printf("THEmx_t::::Print(%s) ==\n",tit);
00111 const double *e = &mHH;
00112 for (int i=0,li=0;i< 5;li+=++i) {
00113 printf("%c ",N[i]);
00114 for (int j=0;j<=i;j++) {
00115 printf("%g\t",e[li+j]);}
00116 printf("\n");
00117 }
00118 }
00119
00120 double THEmx_t::Sign() const
00121 {
00122 const double *E = &mHH;
00123 double dia[5];
00124 for (int i=0,li=0;i<5;li+=++i) {
00125 dia[i] = E[li+i];
00126 if (dia[i]<0) return -(i+1)*11;
00127 for (int j=0;j<=i;j++) {
00128 double dis = dia[i]*dia[j]-E[li+j]*E[li+j];
00129 if (dis<0) return -(i+1+10*(j+1));
00130 } }
00131 return 0;
00132
00133 }
00134
00135
00136 const double Zero = 1.e-6;
00137 static TComplex sgCX1,sgCX2,sgCD1,sgCD2,sgImTet,sgCOne,sgCf1;
00138 static int SqEqu(double *, double *);
00139 #if 0
00140
00141 static int myEqu(double *s, int na, double *b,int nb)
00142 {
00143 StMatrixD mtx(na,na);
00144 double *m = &mtx(1,1);
00145 TCL::trupck(s,m,na);
00146 size_t ierr=0;
00147 mtx.invert(ierr);
00148 if (ierr) return ierr;
00149 for (int ib=0;ib<nb;ib++) {
00150 TCL::vmatl(m,b+ib*na,s,na,na);
00151 memcpy(b+ib*na,s,na*sizeof(*b));
00152 }
00153 TCL::trpck(m,s,na);
00154 return 0;
00155 }
00156
00157 static void Eigen2(const double err[3], double lam[2], double eig[2][2])
00158 {
00159
00160 double spur = err[0]+err[2];
00161 double det = err[0]*err[2]-err[1]*err[1];
00162 double dis = spur*spur-4*det;
00163 if (dis<0) dis = 0;
00164 dis = sqrt(dis);
00165 lam[0] = 0.5*(spur+dis);
00166 lam[1] = 0.5*(spur-dis);
00167 eig[0][0] = 1; eig[0][1]=0;
00168 if (dis>1e-6*spur) {
00169 if (fabs(err[0]-lam[0])>fabs(err[2]-lam[0])) {
00170 eig[0][1] = 1; eig[0][0]= -err[1]/(err[0]-lam[0]);
00171 } else {
00172 eig[0][0] = 1; eig[0][1]= -err[1]/(err[2]-lam[0]);
00173 }
00174 double tmp = sqrt(eig[0][0]*eig[0][0]+eig[0][1]*eig[0][1]);
00175 eig[0][0]/=tmp; eig[0][1]/=tmp;
00176 }
00177 eig[1][0]=-eig[0][1]; eig[1][1]= eig[0][0];
00178 }
00179
00180 static TComplex MyFactor(double rho,double drho,double s)
00181 {
00182
00183
00184
00185
00186
00187 TComplex arr[3],add;
00188 TComplex Sum;
00189 Sum = 0.0;
00190 arr[0] = 1.; arr[1] = 0.;
00191 drho = drho/rho;
00192 double ss = s;
00193 for (int j=2;1;j++) {
00194 arr[2] = -TComplex(0,1)*rho*(arr[1]+drho*arr[0])/double(j);
00195 ss *=s; add = ss*arr[2]; Sum += add;
00196 if (1e-12*Sum.Rho2() > add.Rho2()) break;
00197
00198 arr[0]=arr[1]; arr[1]=arr[2];
00199 }
00200 return Sum;
00201 }
00202 #endif //0
00203
00204 ClassImp(THelixTrack)
00205
00206 THelixTrack::THelixTrack(const double *xyz,const double *dir,double rho
00207 ,double drho)
00208 {
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219
00220
00221
00222
00223
00224
00225
00226
00227
00228
00229
00230
00231 fEmx=0;
00232 Set(xyz,dir,rho,drho);
00233 }
00234
00235 THelixTrack &THelixTrack::operator=(const THelixTrack &from)
00236 {
00237 THEmx_t *save = fEmx;
00238 memcpy(fBeg,from.fBeg,fEnd-fBeg);
00239 fEmx=save;
00240 if (from.fEmx) SetEmx(from.fEmx->Arr());
00241 return *this;
00242 }
00243
00244 THelixTrack::THelixTrack(const THelixTrack &from)
00245 {
00246 fEmx=0;
00247 *this = from;
00248 }
00249
00250 THelixTrack::THelixTrack(const THelixTrack *fr)
00251 {
00252 fEmx=0;
00253 Set(fr->fX,fr->fP,fr->fRho);
00254 }
00255
00256 THelixTrack::~THelixTrack()
00257 { delete fEmx;fEmx=0;}
00258
00259 THelixTrack::THelixTrack()
00260 {
00261 memset(fBeg,0,fEnd-fBeg);
00262 }
00263
00264 void THelixTrack::Set(const double *xyz,const double *dir,double rho
00265 ,double drho)
00266 {
00267 fX[0] = xyz[0]; fX[1] = xyz[1]; fX[2] = xyz[2];
00268 fP[0] = dir[0]; fP[1] = dir[1]; fP[2] = dir[2];
00269 fRho = rho; fDRho=drho;
00270 Build();
00271 }
00272
00273 void THelixTrack::SetEmx(const double* err2xy,const double* err2sz)
00274 {
00275 if (!fEmx) fEmx = new THEmx_t;
00276 fEmx->Set(err2xy,err2sz);
00277 }
00278
00279 void THelixTrack::SetEmx(const double* err)
00280 {
00281 if (!fEmx) fEmx = new THEmx_t;
00282 fEmx->Set(err);
00283 }
00284
00285 void THelixTrack::StiEmx(double err[21]) const
00286 {
00287 enum {kXX
00288 ,kYX,kYY
00289 ,kZX,kZY,kZZ
00290 ,kEX,kEY,kEZ,kEE
00291 ,kPX,kPY,kPZ,kPE,kPP
00292 ,kTX,kTY,kTZ,kTE,kTP,kTT
00293 ,kLN
00294 };
00295 memset(err,0,sizeof(err[0])*kLN);
00296 double cosCA = fP[0]/fCosL;
00297 err[kYY] = fEmx->mHH/(cosCA*cosCA);
00298 err[kZY] = fEmx->mHZ/(cosCA);
00299 err[kZZ] = fEmx->mZZ;
00300 err[kEY] = fEmx->mHA/cosCA;
00301 err[kEZ] = fEmx->mAZ;
00302 err[kEE] = fEmx->mAA;
00303 err[kPY] = fEmx->mHC/cosCA;
00304 err[kPZ] = fEmx->mCZ;
00305 err[kPE] = fEmx->mAC;
00306 err[kPP] = fEmx->mCC;
00307 err[kTY] = fEmx->mHL/(cosCA*fCosL*fCosL);
00308 err[kTZ] = fEmx->mZL/( fCosL*fCosL);
00309 err[kTE] = fEmx->mAL/( fCosL*fCosL);
00310 err[kTP] = fEmx->mCL/( fCosL*fCosL);
00311 err[kTT] = fEmx->mLL/( fCosL*fCosL*fCosL*fCosL);
00312 }
00313
00314 void THelixTrack::Set(double rho,double drho)
00315 {
00316 fRho = rho; fDRho=drho;
00317 }
00318
00319 void THelixTrack::Backward()
00320 {
00321
00322 double d[3];
00323 for (int i=0;i<3;i++) { d[i]=-fP[i];}
00324 Set(fX,d,-fRho,-fDRho);
00325 if(fEmx) fEmx->Backward();
00326 }
00327
00328 void THelixTrack::GetSpot(const double axis[3][3],double emx[3]) const
00329 {
00334
00335
00336 double my[3][3] = {{-fP[1]/fCosL, 0,fP[0]}
00337 ,{ fP[0]/fCosL, 0,fP[1]}
00338 ,{ 0, 1,fP[2]}};
00339
00340 double T[3][3],tmp[3][3],g[6],t[2][2];
00341 TCL::mxmpy (axis[0],my[0],T[0],3,3,3);
00342
00343 TCL::traat(axis[0],g,3,3);
00344 if (fabs(g[0]-1)+fabs(g[1])+fabs(g[2]-1)
00345 +fabs(g[3])+fabs(g[4])+fabs(g[5]-1)>1e-10) {
00346 TCL::trsinv(g,g,3);
00347 memcpy(tmp[0],T[0],sizeof(T));
00348 TCL::trsa (g,tmp[0],T[0],3,3);
00349 }
00350 TCL::vlinco(T[0],1.,T[2],-T[0][2]/T[2][2],t[0],2);
00351 TCL::vlinco(T[1],1.,T[2],-T[1][2]/T[2][2],t[1],2);
00352 double myerr[3]={fEmx->mHH,fEmx->mHZ,fEmx->mZZ};
00353 TCL::trasat(t[0],myerr,emx,2,2);
00354 return;
00355 }
00356
00357 void THelixTrack::Build()
00358 {
00359
00360 double tmp;
00361
00362 tmp = fP[0]*fP[0]+ fP[1]*fP[1]+ fP[2]*fP[2];
00363 if (fabs(tmp-1.) > 1.e-12) {
00364 tmp = ::sqrt(tmp); fP[0] /=tmp; fP[1] /=tmp; fP[2] /=tmp; }
00365
00366 fCosL = ::sqrt(fP[0]*fP[0]+fP[1]*fP[1]);
00367 }
00368
00369 void THelixTrack::MakeMtx(double step,double F[5][5])
00370 {
00371
00372 enum {kH=0,kA,kC,kZ,kL};
00373
00374 double S = step*fCosL;
00375 memset(F[0],0,sizeof(F[0][0])*5*5);
00376
00377 F[kH][kH] = sgCf1.Re()+1.;
00378 double dSdH = sgCf1.Im();
00379
00380 F[kH][kA] = S*sgCOne.Re();
00381 double dSdA = S*sgCOne.Im();
00382
00383 TComplex llCOneD = S*S*expOneD(-sgImTet);
00384 F[kH][kC] = llCOneD.Re();
00385 double dSdC = llCOneD.Im();
00386
00387 F[kA][kH] = -dSdH*fRho;
00388 F[kA][kA] = 1-dSdA*fRho;
00389 F[kA][kC] = S+dSdC*fRho;
00390 F[kC][kC] = 1;
00391
00392 double tanL = fP[2]/fCosL;
00393
00394 F[kZ][kH] = -dSdH*tanL;
00395 F[kZ][kA] = -dSdA*tanL;
00396 F[kZ][kC] = dSdC*tanL;
00397 F[kZ][kZ] = 1;
00398 F[kZ][kL] = S/(fCosL*fCosL);
00399 F[kL][kL] = 1;
00400 }
00401
00402 double THelixTrack::Move(double step)
00403 {
00404 double F[5][5];
00405 Eval(step,fX,fP,fRho);
00406 if (fEmx && fEmx->mHH>0) {
00407 MakeMtx(step,F);
00408 fEmx->Move(F);
00409 }
00410 return step;
00411 }
00412
00413 double THelixTrack::Move(double step,double F[5][5])
00414 {
00415 double xyz[3],dir[3],rho;
00416 Eval(step,xyz,dir,rho);
00417 Set(xyz,dir,rho,fDRho);
00418 MakeMtx(step,F);
00419 if (fEmx && fEmx->mHH>0) fEmx->Move(F);
00420 return step;
00421 }
00422
00423
00424 double THelixTrack::Step(double stmax,const double *surf, int nsurf,
00425 double *xyz, double *dir, int nearest) const
00426 {
00427 int i;
00428 double s[10]={0,0,0,0,0,0,0,0,0,0},tmp=0;
00429 memcpy(s,surf,nsurf*sizeof(surf[0]));
00430
00431 for(i=1;i<nsurf;i++) if (fabs(s[i])>tmp) tmp = fabs(s[i]);
00432 if(fabs(tmp-1.)>0.1) {for(i=0;i<nsurf;i++) s[i]/=tmp;}
00433 double stmin = (nearest)? -stmax:0;
00434
00435
00436
00437 return Step(stmin,stmax,s,nsurf,xyz,dir,nearest);
00438 }
00439
00440
00441
00442 double THelixTrack::Step(double stmin,double stmax, const double *s, int nsurf,
00443 double *xyz, double *dir, int nearest) const
00444 {
00445 int ix,jx,nx,ip,jp;
00446 double poly[4][3],tri[3],sol[2],cos1t,f1,f2,step,ss;
00447 const double *sp[4][4] = {{s+0,s+1,s+2,s+3}, {s+1,s+4,s+7,s+9},
00448 {s+2,s+7,s+5,s+8}, {s+3,s+9,s+8,s+6}};
00449
00450 double myMax = 1./(fabs(fRho*fCosL)+1.e-10);
00451 THelixTrack th(fX,fP,fRho);
00452 cos1t = 0.5*fRho*fCosL;
00453 double totStep=0;
00454 while (2005) {
00455 double hXp[3]={-th.fP[1],th.fP[0],0};
00456 poly[0][0]=1.;poly[0][1]=0.;poly[0][2]=0.;
00457 tri[0]=tri[1]=tri[2]=0;
00458 for(ix=1;ix<4;ix++) {
00459 poly[ix][0] =th.fX [ix-1];
00460 poly[ix][1] =th.fP [ix-1];
00461 poly[ix][2] =hXp[ix-1]*cos1t;
00462 }
00463
00464 nx = (nsurf<=4) ? 1:4;
00465 for(ix=0;ix<nx;ix++) {
00466 for(jx=ix;jx<4;jx++) {
00467 ss = *sp[ix][jx]; if(!ss) continue;
00468 for (ip=0;ip<3;ip++) {
00469 f1 = poly[ix][ip]; if(!f1) continue;
00470 f1 *=ss;
00471 for (jp=0;jp+ip<3;jp++) {
00472 f2 = poly[jx][jp]; if(!f2) continue;
00473 tri[ip+jp] += f1*f2;
00474 } } } }
00475
00476 int nsol = SqEqu(tri,sol);
00477 step = 1.e+12;
00478 if (nsol<0) return step;
00479
00480 if (nearest && nsol>1) {
00481 if(fabs(sol[0])>fabs(sol[1])) sol[0]=sol[1];
00482 nsol = 1;
00483 }
00484 if (nsol) step = sol[0];
00485 if (step < stmin && nsol > 1) step = sol[1];
00486 if (step < stmin || step > stmax) {
00487 nsol = 0;
00488 if (step>0) {step = stmax; stmin+=myMax/2;}
00489 else {step = stmin; stmax-=myMax/2;}}
00490
00491 if (!nsol && fabs(step) < 0.1*myMax) return 1.e+12;
00492 if (fabs(step)>myMax) {step = (step<0)? -myMax:myMax; nsol=0;}
00493
00494 double x[3],d[3];
00495 th.Step(step,x,d);
00496 if (nsol) {
00497 ss = s[0]+s[1]*x[0]+s[2]*x[1]+s[3]*x[2];
00498 if (nsurf > 4) ss += s[4]*x[0]*x[0]+s[5]*x[1]*x[1]+s[6]*x[2]*x[2];
00499 if (nsurf > 7) ss += s[7]*x[0]*x[1]+s[8]*x[1]*x[2]+s[9]*x[2]*x[0];
00500 if (fabs(ss)<1.e-7) {
00501 if (xyz) memcpy(xyz,x,sizeof(*xyz)*3);
00502 if (dir) memcpy(dir,d,sizeof(*dir)*3);
00503 return totStep+step;
00504 } }
00505
00506 stmax -=step; stmin -=step;
00507 if (stmin>=stmax) return 1.e+12;
00508 totStep+=step;
00509 th.Move(step);
00510 }
00511
00512 }
00513
00514
00515 double THelixTrack::StepHZ(const double *su, int nsurf,
00516 double *xyz, double *dir,int nearest) const
00517 {
00518 double tri[3] = {0,0,0};
00519 double f0,fc,fs,R,tet,tet0,tet1,tet2,costet,su45=0,fcs;
00520
00521
00522 R = 1./fRho/fCosL;
00523
00524 f0 = fX[0] - fP[1]*R;
00525 fc = fP[1]*R;
00526 fs = fP[0]*R;
00527
00528 tri[0] = su[0] + su[1]*f0;
00529 tri[1] = su[1]*fc;
00530 tri[2] = su[1]*fs;
00531 if (nsurf >4) {
00532 su45 = 0.5*(su[4]+su[5]);
00533 fcs = fc*fc + fs*fs;
00534 tri[0] += su45*f0*f0 + su45*fcs;
00535 tri[1] += su45*2*f0*fc;
00536 tri[2] += su45*2*f0*fs;
00537 }
00538
00539 f0 = fX[1] + fP[0]*R;
00540 fc = -fP[0]*R;
00541 fs = fP[1]*R;
00542
00543 tri[0] += su[2]*f0;
00544 tri[1] += su[2]*fc;
00545 tri[2] += su[2]*fs;
00546
00547 if (nsurf >4) {
00548 tri[1] += su45*2*f0*fc;
00549 tri[2] += su45*2*f0*fs;
00550 }
00551 costet = -tri[0]/::sqrt(tri[1]*tri[1]+tri[2]*tri[2]);
00552 if(fabs(costet)>1.) return 1.e+12;
00553 tet0 = atan2(tri[2],tri[1]);
00554 tet = acos(costet);
00555 tet1 = tet + tet0;
00556 tet2 = -tet + tet0;
00557
00558 if (tet1 > 2*M_PI) tet1 -= 2*M_PI;
00559 if (tet2 > 2*M_PI) tet2 -= 2*M_PI;
00560 if (nearest) {
00561 if (fabs(tet1)>fabs(tet1-2*M_PI)) tet1 -=2*M_PI;
00562 if (fabs(tet1)>fabs(tet1+2*M_PI)) tet1 +=2*M_PI;
00563 if (fabs(tet2)>fabs(tet2-2*M_PI)) tet2 -=2*M_PI;
00564 if (fabs(tet2)>fabs(tet2+2*M_PI)) tet2 +=2*M_PI;
00565 if (fabs(tet1)>fabs(tet2) ) tet1 =tet2;
00566 return Step(tet1*R,xyz,dir);
00567 } else {
00568 double s1 = tet1*R;
00569 if (s1<=0) s1 += 2*M_PI*fabs(R);
00570 double s2 = tet2*R;
00571 if (s2<=0) s2 += 2*M_PI*fabs(R);
00572 if (s1>s2) s1=s2;
00573 return Step(s1,xyz,dir);
00574 }
00575
00576 }
00577
00578 double THelixTrack::Path(const THelixTrack &th,double *s2) const
00579 {
00580 double SxyMe,SxyHe,SMe,SHe,dSMe=0,dSHe=0;
00581 TCircle tcMe,tcHe;
00582 Fill(tcMe);th.Fill(tcHe);
00583 SxyMe = tcMe.Path(tcHe,&SxyHe);
00584 if (SxyMe>1e33) return 3e33;
00585 THelixTrack thMe(fX,fP,fRho),thHe(th.fX,th.fP,th.fRho);
00586
00587 SMe = SxyMe/thMe.GetCos(); thMe.Move(SMe);
00588 SHe = SxyHe/thHe.GetCos(); thHe.Move(SHe);
00589 for (int ix=0;ix<3;ix++) {
00590 if (fabs(thMe.Pos()[ix]-thHe.Pos()[ix])>100)return 3e33;}
00591 for (int iter=0;iter<10; iter++) {
00592 dSMe = thMe.Path(thHe.Pos());
00593 if (dSMe>1e33) return 3e33;
00594 SMe+=dSMe; thMe.Move(dSMe);
00595 dSHe = thHe.Path(thHe.Pos());
00596 if (dSHe>1e33) return 3e33;
00597 SHe+=dSHe; thHe.Move(dSHe);
00598 if(fabs(dSHe)+fabs(dSMe)<1e-5) break;
00599 }
00600 if(fabs(dSHe)+fabs(dSMe)> 1e-5) return 3e33;
00601 if(s2) *s2=SHe;
00602 return SMe;
00603 }
00604
00605 double THelixTrack::PathX(const THelixTrack &th,double *s2, double *dst, double *xyz) const
00606 {
00607 double ss1,ss2,dd,ss1Best,ss2Best,ddBest=1e33;
00608 double xx[9];
00609 int jkBest=-1;
00610 for (int jk=0;jk<4;jk++) {
00611 THelixTrack th1(this),th2(&th);
00612 if (jk&1) th1.Backward();
00613 if (jk&2) th2.Backward();
00614 ss1 = th1.Path(th2,&ss2);
00615 if (ss1>=1e33) continue;
00616 if (ss2>=1e33) continue;
00617 th1.Eval(ss1,xx+0);
00618 th2.Eval(ss2,xx+3);
00619 TCL::vsub(xx,xx+3,xx+6,3);
00620 dd = TCL::vdot(xx+6,xx+6,3);
00621 if (dd > ddBest) continue;
00622 ddBest = dd; jkBest=jk; ss1Best = ss1; ss2Best = ss2;
00623 if (xyz) TCL::ucopy(xx,xyz,6);
00624 }
00625 if (jkBest<0) { if(s2) *s2=3e33; return 3e33; }
00626 if (jkBest&1) ss1Best = -ss1Best;
00627 if (jkBest&2) ss2Best = -ss2Best;
00628 if (s2 ) *s2 = ss2Best;
00629 if (dst) *dst = ddBest;
00630 return ss1Best;
00631 }
00632
00633 double THelixTrack::Path(double x,double y) const
00634 {
00635 TCircle ht(fX,fP,fRho);
00636 double ar[2]={x,y};
00637 return ht.Path(ar)/fCosL;
00638 }
00639
00640 double THelixTrack::Step(const double *point,double *xyz, double *dir) const
00641 {
00642
00643 static int nCount=0; nCount++;
00644 TComplex cpnt(point[0]-fX[0],point[1]-fX[1]);
00645 TComplex cdir(fP[0],fP[1]); cdir /=TComplex::Abs(cdir);
00646 double step[3]={0,0,0};
00647
00648
00649 int zStep=0;
00650 if (fabs(fP[2]) > 0.01){
00651 zStep = 1;
00652 step[1] = (point[2]-fX[2])/fP[2];
00653 }
00654
00655
00656 {
00657 cpnt /= cdir;
00658 if (fabs(cpnt.Re()*fRho) < 0.01) {
00659 step[2]=cpnt.Re();
00660 } else {
00661 double rho = fRho;
00662 for (int i=0;i<2;i++) {
00663 TComplex ctst = (1.+TComplex(0,1)*rho*cpnt);
00664 ctst /=TComplex::Abs(ctst);
00665 ctst = TComplex::Log(ctst);
00666 step[2]= ctst.Im()/rho;
00667 if (!fDRho) break;
00668 rho = fRho+ 0.5*fDRho*step[2];
00669 }
00670 }
00671 step[2]/=fCosL;
00672 }
00673
00674 if (zStep) {
00675 double p = GetPeriod();
00676 int nperd = (int)((step[1]-step[2])/p);
00677 if (step[2]+nperd*p>step[1]) nperd--;
00678 if (fabs(step[2]-step[1]+(nperd+0)*p)
00679 >fabs(step[2]-step[1]+(nperd+1)*p)) nperd++;
00680 step[2]+=(nperd)*p;
00681 }
00682 step[0] = step[2];
00683
00684 double ds = step[1]-step[2];
00685 if (zStep && fabs(ds)>1.e-5) {
00686 double dz = ds*fP[2];
00687 step[0] += dz*dz/ds;
00688 }
00689
00690
00691 double xnear[6],ss=0; double* pnear=xnear+3;
00692
00693 double dstep = 1.e+10,oldStep=dstep,dztep;
00694 double lMax = step[0]+0.25*GetPeriod();
00695 double lMin = step[0]-0.25*GetPeriod();
00696
00697 if (zStep) {
00698 lMax = (step[1]>step[2])? step[1]:step[2];
00699 lMin = (step[1]>step[2])? step[2]:step[1];}
00700 int iter=99,icut=1;
00701 THelixTrack local(this);
00702 local.Move(step[0]);
00703 lMax-=step[0];lMin-=step[0];
00704 local.Step(0.,xnear,pnear);
00705 for (; iter; iter--)
00706 {
00707 double diff = (icut)? lMax-lMin: fabs(dstep);
00708 if (diff < 0.1) {
00709 if (diff < 1.e-6) break;
00710 double tmpxy = fabs(point[0]-xnear[0])+fabs(point[1]-xnear[1]);
00711 double tmpz = fabs(point[2]-xnear[2]);
00712 if (fabs(fCosL) *diff <tmpxy*1e-6
00713 &&fabs(pnear[2])*diff <tmpz *1e-6) break;
00714 if (tmpxy+tmpz < 1.e-6) break;
00715 }
00716
00717 local.Step(ss,xnear,pnear);
00718 dstep = 0; icut = 0;
00719 for (int i=0;i<3;i++) {dstep += pnear[i]*(point[i]-xnear[i]);}
00720 if(dstep<0) {
00721 lMax = ss; dztep = -0.5*(lMax-lMin);
00722 if (dstep<dztep || fabs(dstep)>0.7*oldStep) {icut=1;dstep = dztep;}
00723 } else {
00724 lMin = ss; dztep = 0.5*(lMax-lMin);
00725 if (dstep>dztep || fabs(dstep)>0.7*oldStep) {icut=1;dstep = dztep;}
00726 }
00727 ss += dstep;
00728 oldStep=fabs(dstep);
00729 }
00730
00731 if (!iter){ printf("*** Problem in THElixTrack::Step(vtx) ***\n");
00732 printf("double vtx[3]={%g,%g,%g};",point[0],point[1],point[2]);
00733 Print();}
00734 assert(iter);
00735 step[0]+=ss;
00736 return (xyz) ? Step(step[0],xyz,dir) : step[0];
00737 }
00738
00739 double THelixTrack::Dca(const double *point,double *dcaErr) const
00740 {
00741 double x[3],T[3][3],emx[3];
00742 double s = Path(point,x,T[2]);
00743 for (int i=0;i<3;i++) {T[0][i]=point[i]-x[i];}
00744 double dca = sqrt(T[0][0]*T[0][0]+T[0][1]*T[0][1]+T[0][2]*T[0][2]);
00745 if (!dcaErr) return dca;
00746
00747 for (int i=0;i<3;i++) {T[0][i]/=dca;}
00748 T[1][0]=T[0][1]*T[2][2]-T[2][1]*T[0][2];
00749 T[1][1]=T[0][2]*T[2][0]-T[2][2]*T[0][0];
00750 T[1][2]=T[0][0]*T[2][1]-T[2][0]*T[0][1];
00751
00752 THelixTrack th(*this);
00753 th.Move(s);
00754 th.GetSpot(T,emx);
00755 *dcaErr=emx[0];
00756 return dca;
00757 }
00758
00759 double THelixTrack::Dca(double x,double y,double *dcaErr) const
00760 {
00761 double dir[3]={fP[0],fP[1],0};
00762 THelixTrack hlx(fX,dir,fRho);
00763 if (fEmx) hlx.SetEmx(fEmx->Arr());
00764 double vtx[3]={x,y,fX[2]};
00765 return hlx.Dca(vtx,dcaErr);
00766 }
00767
00768
00769
00770 double THelixTrack::Dca(const double point[3]
00771 ,double &dcaXY,double &dcaZ,double dcaEmx[3],int kind) const
00779 {
00780 double dif[3];
00781 double s = 0;
00782 assert(kind==2 || kind==3);
00783 if (kind==3) s = Path(point);
00784 else s = Path(point[0],point[1]);
00785
00786 THelixTrack th(*this);
00787 th.Move(s);
00788 const double *x=th.Pos();
00789 const double *d=th.Dir();
00790
00791 for (int i=0;i<3;i++) {dif[i]=x[i]-point[i];}
00792 double nor = th.GetCos();
00793 double T[3][3]={{-d[1]/nor, d[0]/nor, 0}
00794 ,{ 0, 0, 1}
00795 ,{ d[0]/nor, d[1]/nor, 0}};
00796
00797 dcaXY = T[0][0]*dif[0]+T[0][1]*dif[1];
00798 dcaZ = dif[2];
00799 if (!dcaEmx) return s;
00800 THEmx_t *emx =th.Emx();
00801 dcaEmx[0] = emx->mHH;
00802 dcaEmx[1] = 0;
00803
00804 dcaEmx[2] = emx->mZZ*pow(th.GetCos(),4);
00805 return s;
00806 }
00807
00808
00809
00810 double THelixTrack::GetPeriod() const
00811 {
00812 double per = (fabs(fRho) > 1.e-10) ? fabs(2.*M_PI/fRho):1.e+10;
00813 return per/fCosL;
00814 }
00815
00816 void THelixTrack::Rot(double angle)
00817 {
00818 Rot(cos(angle),sin(angle));
00819 }
00820
00821 void THelixTrack::Rot(double cosa,double sina)
00822 {
00823 TComplex CX(fX[0],fX[1]),CP(fP[0],fP[1]);
00824 TComplex A (cosa,sina);
00825 CX *=A; fX[0] = CX.Re(); fX[1]=CX.Im();
00826 CP *=A; fP[0] = CP.Re(); fP[1]=CP.Im();
00827 }
00828
00829 void THelixTrack::Streamer(TBuffer &){}
00830
00831 void THelixTrack::Print(Option_t *) const
00832 {
00833 printf("\n THelixTrack::this = %p\n",(void*)this);
00834 printf(" THelixTrack::fX[3] = { %f , %f ,%f }\n",fX[0],fX[1],fX[2]);
00835 printf(" THelixTrack::fP[3] = { %f , %f ,%f }\n",fP[0],fP[1],fP[2]);
00836 printf(" THelixTrack::fRho = %f \n\n",fRho);
00837
00838 printf("double xyz[3] = {%g,%g,%g};\n" ,fX[0],fX[1],fX[2]);
00839 printf("double dir[3] = {%g,%g,%g};\n" ,fP[0],fP[1],fP[2]);
00840 printf("double Rho = %g;\n" ,fRho);
00841 printf("THelixTrack *ht = new THelixTrack(xyz,dir,Rho);\n");
00842
00843 }
00844
00845 void THelixTrack::TestMtx()
00846 {
00847 enum {kH=0,kA,kC,kZ,kL};
00848 const static char* T="HACZL";
00849 double Dir[4][3],X[4][3]={{0}},Rho[2],step,F[5][5],Del,Dif;
00850 double maxEps = 0;
00851 int nErr=0;
00852 int iR = 10+ gRandom->Rndm()*100;
00853 int iAlf=30+ gRandom->Rndm()*100;
00854 int iLam=10+ gRandom->Rndm()*60;
00855 step = gRandom->Rndm()*6*iR;
00856
00857 Rho[0] = 1./iR;
00858 double alf = iAlf/180.*M_PI;
00859 double lam = iLam/180.*M_PI;
00860 Dir[0][0] = cos(lam)*cos(iAlf/180.*M_PI);
00861 Dir[0][1] = cos(lam)*sin(iAlf/180.*M_PI);
00862 Dir[0][2] = sin(lam);
00863 THelixTrack tc(X[0],Dir[0],Rho[0]);
00864 tc.Eval(step,X[1],Dir[1]);
00865
00866 tc.Move(step,F);
00867
00868 printf("TestMtx: Angle=%d Lam=%d \tRad=%d Step=%d \n",iAlf,iLam,iR,int(step));
00869
00870 for (int iHAR=0;iHAR<5;iHAR++) {
00871 memcpy(X[2] ,X[0] ,sizeof(X[0][0]) *3);
00872 memcpy(Dir[2],Dir[0],sizeof(Dir[0][0])*3);
00873 Del = 0;
00874 Rho[1]=Rho[0];
00875 switch (iHAR) {
00876 case kH: {
00877 Del = 0.001*iR;
00878 X[2][0] += -Dir[0][1]*Del/cos(lam);
00879 X[2][1] += Dir[0][0]*Del/cos(lam);
00880 break;}
00881
00882 case kA: {
00883 Del = M_PI/180*0.01;
00884 Dir[2][0] = cos(lam)*cos(alf+Del);
00885 Dir[2][1] = cos(lam)*sin(alf+Del);
00886 Dir[2][2] = sin(lam);
00887 break;}
00888
00889 case kC: {
00890 Del = Rho[0]*0.005;
00891 Rho[1] = Rho[0]+Del;
00892 break;}
00893 case kZ: {
00894 Del = 0.02;
00895 X[2][2] += Del;
00896 break;}
00897 case kL: {
00898 Del = M_PI/180*0.1;
00899 Dir[2][0] = cos(lam+Del)*cos(alf);
00900 Dir[2][1] = cos(lam+Del)*sin(alf);
00901 Dir[2][2] = sin(lam+Del);
00902 break;}
00903 }
00904
00905 THelixTrack tcc(X[2],Dir[2],Rho[1]);
00906 tcc.Move(step);
00907 double myStep = tcc.Path(X[1][0],X[1][1]);
00908 tcc.Eval(myStep,X[3],Dir[3]);
00909
00910 for (int jHAR=0;jHAR<5; jHAR++) {
00911 if (jHAR==kC) continue;
00912 if (jHAR==kL) continue;
00913 switch(jHAR) {
00914 case kH: {
00915 Dif = (X[3][0]-X[1][0])*(-Dir[1][1])
00916 + (X[3][1]-X[1][1])*( Dir[1][0]);
00917 Dif/=cos(lam);
00918 break;}
00919 case kA: {
00920 Dif = atan2(Dir[3][1],Dir[3][0])
00921 -atan2(Dir[1][1],Dir[1][0]);
00922 if (Dif> M_PI) Dif-=2*M_PI;
00923 if (Dif< -M_PI) Dif+=2*M_PI;
00924 break;}
00925 case kZ: {
00926 Dif = X[3][2]-X[1][2];
00927 break;}
00928 }
00929 double est = Dif/Del;
00930 double eps = fabs(est-F[jHAR][iHAR])*2
00931 /(fabs(est)+fabs(F[jHAR][iHAR]+1e-6));
00932 if (eps>maxEps) maxEps=eps;
00933 if (eps < 1e-2) continue;
00934 nErr++;
00935 printf(" m%c%c \t%g \t%g \t%g\n",
00936 T[jHAR],T[iHAR],F[jHAR][iHAR],est,eps);
00937
00938 } }
00939 printf("TestMtx: %d errors maxEps=%g\n",nErr,maxEps);
00940
00941 }
00942
00943
00944 int SqEqu(double *cba, double *sol)
00945 {
00946
00947
00948
00949
00950
00951
00952
00953
00954
00955
00956
00957
00958
00959
00960
00961
00962
00963
00964 const double zero2=1.e-12;
00965 double swap,a,b,c,amx,dis,bdis;
00966 int nsol;
00967
00968
00969 a = cba[2]; b = cba[1]*0.5; c = cba[0];
00970 if (b < 0.) {a = -a; b = -b; c = -c;}
00971 amx = fabs(a); if (amx<b) amx = b; if (amx<fabs(c)) amx = fabs(c);
00972 if (amx <= 0.) return -1;
00973 a = a/amx; b = b/amx; c = c/amx;
00974
00975 dis = b*b - a*c;
00976 nsol = 1;
00977 if (fabs(dis) <= zero2) dis = 0;
00978 if (dis < 0.) { nsol = 0; dis = 0.;}
00979
00980 dis = ::sqrt(dis); bdis = b + dis;
00981 if (fabs(c) > 1.e+10*bdis) return -1;
00982 sol[0] = 0.;
00983 if (fabs(bdis) <= 0.) return nsol;
00984 sol[0] = (-c/bdis);
00985 if (dis <= 0.) return nsol;
00986 if (bdis >= 1.e+10*fabs(a)) return nsol;
00987 nsol = 2; sol[1] = (-bdis/a);
00988 if (sol[0] > sol[1]) { swap = sol[0]; sol[0] = sol[1]; sol[1] = swap;}
00989 return nsol;
00990 }
00991
00992 double THelixTrack::Eval(double step, double *xyz, double *dir,double &rho) const
00993 {
00994 Eval(step,xyz,dir);
00995 rho = fRho +(step*fCosL)*fDRho;
00996 return step;
00997 }
00998
00999 double THelixTrack::Eval(double step, double *xyz, double *dir) const
01000 {
01001 if (!step) {
01002 if (xyz) memcpy(xyz,fX,sizeof(fX));
01003 if (dir) memcpy(dir,fP,sizeof(fP));
01004 return 0.;
01005 }
01006
01007 double ztep = step*fCosL;
01008 double teta = ztep*(fRho+0.5*ztep*fDRho);
01009
01010 sgCX1 = TComplex(fX[0] ,fX[1]);
01011 sgCD1 = TComplex(fP[0],fP[1])/fCosL;
01012 sgImTet = TComplex(0,teta);
01013 sgCOne = expOne(sgImTet);
01014 sgCf1 = sgImTet*sgCOne;
01015 sgCD2 = sgCD1*sgCf1+sgCD1;
01016 sgCX2 = sgCD1*sgCOne*ztep;
01017
01018 if (xyz) {
01019 xyz[0] = sgCX2.Re()+sgCX1.Re();
01020 xyz[1] = sgCX2.Im()+sgCX1.Im();
01021 xyz[2] = fX[2]+fP[2]*step;
01022 }
01023 if (dir) {
01024 sgCD2/= TComplex::Abs(sgCD2);
01025 dir[0] = sgCD2.Re()*fCosL;
01026 dir[1] = sgCD2.Im()*fCosL;
01027 dir[2] = fP[2];
01028 }
01029 return step;
01030 }
01031
01032 void THelixTrack::Fill(TCircle &circ) const
01033 {
01034 circ.fX[0]=fX[0];
01035 circ.fX[1]=fX[1];
01036 circ.fD[0]=fP[0]/fCosL;
01037 circ.fD[1]=fP[1]/fCosL;
01038 circ.fRho=fRho;
01039 if (fEmx) circ.SetEmx(fEmx->Arr());
01040 }
01041
01042 void THelixTrack::InvertMtx(double F[5][5])
01043 {
01044 static const int minus[][2] = {{0,1},{1,0},{1,2},{3,0},{3,2},{3,4},{-1,0}};
01045 for (int i=0;minus[i][0]>=0;i++) {
01046 F[minus[i][0]][minus[i][1]] = -F[minus[i][0]][minus[i][1]];
01047 }
01048 }
01049
01050
01051
01052
01053
01054
01055
01056
01057
01058
01059
01060
01061 void THelixTrack::Show(double len, const THelixTrack *other) const
01062 {
01063 static TCanvas *myCanvas = 0;
01064 int kolor[2]={kRed,kBlue};
01065
01066 TGraph *ciGraph[2][2] = {{0}};
01067 TGraph *ptGraph[2][2] = {{0}};
01068 TGraph *szGraph[2] = {0};
01069
01070 double x[100],y[100],z[100],l[100],xyz[3];
01071 double X[4],Y[4],Z[4],L[4];
01072 const THelixTrack *th[]={this,other};
01073 int nH = (other)? 2:1;
01074 for (int ih=0;ih<nH;ih++) {
01075 double rho = fabs(th[ih]->GetRho());
01076 double step = 0.01*(1./(rho+1e-10));
01077
01078 if (step>fabs(len)*0.10) step=fabs(len)*0.1;
01079 if (step<fabs(len)*0.01) step=fabs(len)*0.01;
01080
01081
01082 int nPts = (int)(fabs(len)/step);
01083 step = fabs(len)/nPts;
01084 if (len<0) {len = -len; step = -step;}
01085 for (int ipt=0; ipt<nPts; ipt++) {
01086 double s = ipt*step;
01087 th[ih]->Eval(s,xyz);
01088 l[ipt]=s; x[ipt]=xyz[0]; y[ipt]=xyz[1], z[ipt]=xyz[2];
01089 }
01090 ciGraph[ih][0] = new TGraph(nPts , x, y);
01091 ciGraph[ih][1] = new TGraph(nPts , l, z);
01092 ciGraph[ih][0]->SetLineColor(kolor[ih]);
01093 ciGraph[ih][1]->SetLineColor(kolor[ih]);
01094 ptGraph[ih][0] = new TGraph( 1 , x, y);
01095 ptGraph[ih][1] = new TGraph( 1 , l, z);
01096 ptGraph[ih][0]->SetMarkerColor(kolor[ih]);
01097 ptGraph[ih][1]->SetMarkerColor(kolor[ih]);
01098
01099 X[ih*2+0]=x[0]; X[ih*2+1]=x[nPts-1];
01100 Y[ih*2+0]=y[0]; Y[ih*2+1]=y[nPts-1];
01101 Z[ih*2+0]=z[0]; Z[ih*2+1]=z[nPts-1];
01102 L[ih*2+0]=l[0]; L[ih*2+1]=l[nPts-1];
01103 }
01104
01105 szGraph[0] = new TGraph(nH*2 , X, Y);
01106 szGraph[1] = new TGraph(nH*2 , L, Z);
01107
01108 myCanvas = new TCanvas("THelixTrack_Show","",600,800);
01109 myCanvas->Divide(1,2);
01110 for (int ipad=0;ipad<2;ipad++) {
01111 myCanvas->cd(ipad+1);
01112 szGraph[ipad]->Draw("AP");
01113 for (int ih = 0;ih<nH;ih++) {
01114 ptGraph[ih][ipad]->Draw("same *");
01115 ciGraph[ih][ipad]->Draw("same CP");
01116 } }
01117
01118
01119 myCanvas->Modified();
01120 myCanvas->Update();
01121 while(!gSystem->ProcessEvents()){gSystem->Sleep(200);};
01122
01123 }
01124
01125 void THelixTrack::Test1()
01126 {
01127 double surf[4] = {-11.32212856152224, 0.50109792630239824, -0.86539108263698283, 0.00078459561521909921};
01128 double xyz[3] = {-0.0206564,-0.0153429,0.0285582};
01129 double dir[3] = {0.450295,-0.596426,-0.664463};
01130 double Rho = 0.00678696;
01131 THelixTrack TH(xyz,dir,Rho);
01132
01133 double s = TH.Step(100000,surf,4);
01134 printf("s=%g = 15.3589\n",s);
01135 }
01136
01137 void THelixTrack::Test2()
01138 {
01139 double diff[3];
01140
01141 double xyz[3] = {-60.0301,1.51445,-1.57283};
01142 double dir[3] = {-0.849461,0.526419,0.0360391};
01143 double Rho = 0.00363571;
01144 THelixTrack ht(xyz,dir,Rho);
01145
01146 double MyHit[3]= { -177.673, 41.305, 2.90798};
01147 double MyClo[3];
01148
01149 printf("%s= %g %g %g\n","MyHit",MyHit[0],MyHit[1],MyHit[2]);
01150 double s = ht.Step(MyHit,MyClo);
01151 ht.Step(s,MyClo);
01152 TCL::vsub(MyClo,MyHit,diff,3);
01153 double MyDist = sqrt(TCL::vdot(diff,diff,3));
01154 printf("%s= %g %g %g\n","MyClo ",MyClo[0],MyClo[1],MyClo[2]);
01155 printf("MustBe= -177.661 41.4145 2.94559\n");
01156
01157 printf("%s= %g %g %g\n","MyDif ",diff[0],diff[1],diff[2]);
01158 printf("MustBe= 0.0122709 0.109539 0.0376077\n");
01159 printf("%s=%g\n","MyS ",s);
01160 printf("MustBe=125.375\n");
01161 printf("%s= %g\n","MyDist",MyDist);
01162 printf("MustBe= 0.116463\n");
01163 }
01164
01165 void THelixTrack::Test3()
01166 {
01167 double xyz[3] = {100,200,300};
01168 double dir[3] = {-0.224845,-0.491295,-0.841471};
01169 double Rho = 0.02;
01170 double sur[8]={-120,1,0,0,0,0,0};
01171 THelixTrack *ht = new THelixTrack(xyz,dir,Rho);
01172 double newX[3],newD[3];
01173 ht->Backward();
01174 double s = ht->Step(1000.,sur,4,newX,newD);
01175 printf("Result: s=%g newX=(%g %g %g) newD=(%g %g %g)\n"
01176 ,s,newX[0],newX[1],newX[2],newD[0],newD[1],newD[2]);
01177
01178 printf("MustBe: s=56.1931 newX=(120 222.222 347.285) newD=(0.464979 0.275174 0.841471)\n\n");
01179
01180 sur[6]=1e-6;
01181 s = ht->Step(1000.,sur,7,newX,newD);
01182 printf("Result: s=%g newX=(%g %g %g) newD=(%g %g %g)\n"
01183 ,s,newX[0],newX[1],newX[2],newD[0],newD[1],newD[2]);
01184 printf("MustBe: s=55.9338 newX=(119.88 222.151 347.067) newD=(0.464206 0.276476 0.841471)\n\n");
01185 }
01186
01187 void THelixTrack::Test4()
01188 {
01189 double surf[7] = {-100, 0, 0, 0, 1,1,0};
01190 double xyz[3] = {-0.0206564,-0.0153429,0.0285582};
01191 double dir[3] = {0.450295,-0.596426,-0.664463};
01192 double Rho = 0.00678696;
01193 double xnew[3];
01194 THelixTrack TH(xyz,dir,Rho);
01195
01196 double s = TH.Step(100000,surf,7,xnew);
01197 double dif = xnew[0]*xnew[0]+xnew[1]*xnew[1]-100;
01198 printf("s=%g dif=%g\n",s,dif);
01199
01200 }
01201
01202 void THelixTrack::Test5()
01203 {
01204 double pars[4][2][7] = {
01205 {{0,0,0, 1,1,1, -0.001},{0,0, 0, -1,1,-1,0.002}},
01206 {{0,0,1, 1,1,1, -0.001},{0,0,-1, -1,1,-1,0.002}},
01207 {{0,0,1, 1,1,1, -0.001},{0,0,-1, 1,1,-1,0.002}},
01208 };
01209 for (int ip=0;ip<3;ip++) {
01210 THelixTrack th1(pars[ip][0],pars[ip][0]+3,pars[ip][0][6]);
01211 THelixTrack th2(pars[ip][1],pars[ip][1]+3,pars[ip][1][8]);
01212 th1.Move(-50);
01213 th2.Move(-50);
01214 double s2;
01215 double s1 = th1.Path(th2,&s2);
01216 th1.Move(s1);
01217 th2.Move(s2);
01218 double diff[3];
01219 TCL::vsub(th1.Pos(),th2.Pos(),diff,3);
01220 double dist = sqrt(TCL::vdot(diff,diff,3));
01221 printf("s1=%g s2=%g dist=%g\n",s1,s2,dist);
01222
01223 }
01224 }
01225
01226 void THelixTrack::TestDer()
01227 {
01228 enum {kH=0,kA,kC,kZ,kL};
01229 double fak = 0.1;
01230 double D[5][3],X[5][3]={{0}},Rho[5],step,F[5][5];
01231 int iR = 10+ gRandom->Rndm()*100;
01232 int iAlf=30+ gRandom->Rndm()*100;
01233 int iLam=10+ gRandom->Rndm()*60;
01234 step = gRandom->Rndm()*6*iR;
01235
01236 Rho[0] = 1./iR;
01237 double alf = iAlf/180.*M_PI;
01238 double lam = iLam/180.*M_PI;
01239 D[0][0] = cos(lam)*cos(alf);
01240 D[0][1] = cos(lam)*sin(alf);
01241 D[0][2] = sin(lam);
01242 THelixTrack hlx0(X[0],D[0],Rho[0]);
01243 THelixTrack hlx1(hlx0);
01244 hlx1.Move(step,F);
01245 hlx1.Eval(0,X[2],D[2]);
01246 Rho[2]=hlx1.GetRho();
01247
01248
01249 printf("TestDer: Angle=%d Lam=%d \tRad=%d Step=%d \n",iAlf,iLam,iR,int(step));
01250 hlx0.Eval(0,X[1],D[1]); Rho[1]=hlx0.GetRho();
01251
01252 double dH = iR *0.01 *(gRandom->Rndm()-0.5)*fak;
01253 double dAlf = M_PI/180*0.1 *(gRandom->Rndm()-0.5)*fak;
01254 double dRho = Rho[0] *0.1 *(gRandom->Rndm()-0.5)*fak;
01255 double dZ = (gRandom->Rndm()-0.5)*fak;
01256 double dLam = M_PI/180*0.1 *(gRandom->Rndm()-0.5)*fak;
01257 double dA[5] = {dH,dAlf,dRho,dZ ,dLam}, dB[5];
01258 D[1][0] = cos(lam+dLam)*cos(alf+dAlf);
01259 D[1][1] = cos(lam+dLam)*sin(alf+dAlf);
01260 D[1][2] = sin(lam+dLam);
01261 X[1][0] += -D[0][1]*dH/cos(lam);
01262 X[1][1] += D[0][0]*dH/cos(lam);
01263 X[1][2] += dZ;
01264 Rho[1]+=dRho;
01265 THelixTrack hlxM(X[1],D[1],Rho[1]);
01266 hlxM.Move(step);
01267 double dL = hlxM.Path(X[2][0],X[2][1]);
01268 hlxM.Move(dL);
01269 hlxM.Eval(0,X[3],D[3]); Rho[3]=hlxM.GetRho();
01270
01271 TCL::vmatl(F[0],dA,dB,5,5);
01272
01273
01274 memcpy(X[4],X[2],3*sizeof(double));
01275 memcpy(D[4],D[2],3*sizeof(double));
01276 Rho[4]=Rho[2];
01277 double myAlf = atan2(D[2][1],D[2][0]);
01278 double myLam = asin(D[2][2]) ;
01279 D[4][0] = cos(myLam+dB[kL])*cos(myAlf+dB[kA]);
01280 D[4][1] = cos(myLam+dB[kL])*sin(myAlf+dB[kA]);
01281 D[4][2] = sin(myLam+dB[kL]);
01282 X[4][0] += -D[2][1]*dB[kH]/cos(myLam);
01283 X[4][1] += D[2][0]*dB[kH]/cos(myLam);
01284 X[4][2] += dB[kZ];
01285 Rho[4] += dB[kC];
01286
01287 THelixTrack hlxD(X[4],D[4],Rho[4]);
01288 dL = hlxD.Path(X[2][0],X[2][1]);
01289 hlxD.Move(dL);
01290 hlxD.Eval(0,X[4],D[4]);
01291 Rho[4]=hlxD.GetRho();
01292
01293
01294
01295
01296 double dC[5];
01297 dC[kH] = (-(X[3][0]-X[2][0])*D[2][1]+(X[3][1]-X[2][1])*D[2][0])/cos(myLam);
01298 dC[kA] = atan2(D[3][1],D[3][0])-atan2(D[2][1],D[2][0]);
01299 dC[kC] = Rho[3]-Rho[2];
01300 dC[kZ] = X[3][2]-X[2][2];
01301 dC[kL] = asin(D[3][2])-asin(D[2][2]);
01302 for (int i=0;i<5;i++) {printf(" %d - %g == %g\n",i,dB[i],dC[i]);}
01303
01304 double myDelta = 0;
01305 for (int i=0;i<3;i++){myDelta+=pow(X[4][i]-X[3][i],2);}
01306 myDelta = sqrt(myDelta);
01307
01308 double ihDelta = 0;
01309 for (int i=0;i<3;i++){ihDelta+=pow(X[2][i]-X[3][i],2);}
01310 ihDelta = sqrt(ihDelta);
01311 printf( "\n*** DELTA = %g << %g ***\n",myDelta,ihDelta);
01312
01313 }
01314
01315
01316
01317
01318
01319 ClassImp(TCircle)
01320
01321 TCircle::TCircle()
01322 {
01323 Set(0,0,0.);
01324 fEmx = 0;
01325 }
01326
01327 TCircle::~TCircle()
01328 {delete fEmx;fEmx=0;}
01329
01330
01331 TCircle::TCircle(const double *x,const double *d,double rho)
01332 {
01333 Set(x,d,rho);
01334 fEmx = 0;
01335 }
01336
01337 void TCircle::Set(const double *x,const double *d,double rho)
01338 {
01339 fX[0]=0; fX[1]=0; fD[0]=0; fD[1]=0;
01340 if (x) {fX[0]=x[0];fX[1]=x[1];}
01341 if (d) {
01342 fD[0]=d[0];fD[1]=d[1];
01343 double n = sqrt(fD[0]*fD[0]+fD[1]*fD[1]);
01344 fD[0]/=n; fD[1]/=n;
01345 }
01346 fRho = rho;
01347 }
01348
01349 TCircle::TCircle(const TCircle& fr)
01350 {
01351 fEmx = 0;
01352 *this = fr;
01353 }
01354
01355 TCircle::TCircle(const TCircle* fr)
01356 {
01357 fEmx = 0;
01358 Set(fr->fX,fr->fD,fr->fRho);
01359 }
01360
01361 TCircle &TCircle::operator=(const TCircle& fr)
01362 {
01363 Set(fr.fX,fr.fD,fr.fRho);
01364 if (fr.fEmx) SetEmx(fr.fEmx->Arr());
01365 return *this;
01366 }
01367
01368
01369 void TCircle::Clear(const char *)
01370 {
01371 if (fEmx) fEmx->Clear();
01372 memset(fX,0,(char*)(&fEmx)-(char*)fX);
01373 }
01374
01375
01376
01377 void TCircle::SetEmx(const double *err)
01378 {
01379 if(!fEmx) fEmx = new TCEmx_t;
01380 fEmx->Set(err);
01381 }
01382
01383 void TCircle::Nor(double *norVec) const
01384 {
01385
01386
01387 norVec[0] = fD[1]; norVec[1] = -fD[0];
01388 if (fRho>=0) return;
01389 norVec[0] = -norVec[0];norVec[1] = -norVec[1];
01390 }
01391
01392 void TCircle::GetCenter(double cent[3]) const
01393 {
01394 double R;
01395 if (fabs(fRho) > 1e-10) { R = 1./fRho ;}
01396 else { R =( fRho>0) ? 1e10:-1e10;}
01397 cent[0] = fX[0]-fD[1]*R;
01398 cent[1] = fX[1]+fD[0]*R;
01399 }
01400
01401 void TCircle::Print(const char* txt) const
01402 {
01403 if (!txt) txt="";
01404 printf("TCircle(%s): x,y=%g %g dir=%g %g curv=%g\n",txt,fX[0],fX[1],fD[0],fD[1],fRho);
01405 if (!fEmx) return;
01406 printf("Errs: %g\n" ,fEmx->mHH);
01407 printf(" : %g \t%g\n" ,fEmx->mHA,fEmx->mAA);
01408 printf(" : %g \t%g \t%g\n",fEmx->mHC,fEmx->mAC,fEmx->mCC);
01409 }
01410
01411 double TCircle::Path(const double pnt[2]) const
01412 {
01413 TComplex CX1(pnt[0]-fX[0],pnt[1]-fX[1]);
01414 TComplex CP(fD[0],fD[1]);
01415 TComplex CXP = TComplex(0,1)*CX1/CP;
01416 TComplex CXPRho = CXP*fRho;
01417 double s;
01418 if (TComplex::Abs(CXPRho)>0.001) {
01419 s = TComplex::Log(1.+CXPRho).Im()/fRho;
01420 } else {
01421 s = (CXP*(1.-CXPRho*(0.5-CXPRho*(1/3.-CXPRho*0.25)))).Im();
01422 }
01423
01424
01425
01426
01427 return s;
01428 }
01429
01430 double TCircle::Path(const double *pnt,const double *exy) const
01431 {
01432 double x[2],d[2];
01433 double s = Path(pnt);
01434 if (!exy || exy[0]<=0) return s;
01435 Eval(s,x,d);
01436 double k = (x[0]-pnt[0])*(d[1]) + (x[1]-pnt[1])*(-d[0]);
01437 double t =((d[1]*d[1]-d[0]*d[0])*exy[1]-d[1]*d[0]*(exy[2]-exy[0]))
01438 /( d[0]*d[0]*exy[2] -2*d[1]*d[0]*exy[1]+d[1]*d[1]*exy[0]);
01439
01440 return s+k*t;
01441 }
01442
01443 double TCircle::Path(const TCircle &hlx,double *S2) const
01444 {
01445 double s,s1=3e33,s2=3e33;
01446 const static TComplex Im(0.,1.);
01447 const TCircle *one = this;
01448 const TCircle *two = &hlx;
01449 if (fabs(fRho) > fabs(hlx.fRho)) { one=two; two=this; }
01450 double Rho1 = one->Rho();
01451 double Rho2 = two->Rho();
01452 int kase = 0;
01453 if (fabs(Rho1) > 1e-4) kase+=1;
01454 if (fabs(Rho2) > 1e-4) kase+=2;
01455
01456 int nSol = 0;
01457 TComplex CX1(one->Pos()[0],one->Pos()[1]);
01458 TComplex CX2(two->Pos()[0],two->Pos()[1]);
01459 TComplex CP1(one->Dir()[0],one->Dir()[1]);
01460 TComplex CP2(two->Dir()[0],two->Dir()[1]);
01461 TComplex CdX = CX2-CX1;
01462 double L[2];
01463 switch(kase) {
01464 case 2:;
01465 case 3: {
01466 if (kase==3) {
01467 TComplex A = CP1/CP2*(Rho2/Rho1);
01468 TComplex B = 1.-CdX/CP2*(Im*Rho2) - CP1/CP2*(Rho2/Rho1);
01469 double a = A.Rho();
01470 double b = B.Rho();
01471 double alfa = A.Theta();
01472 double beta = B.Theta();
01473 double myCos = (1.-(a*a+b*b))/(2.*a*b);
01474 double myAng = 0;
01475 if (myCos>= 1.) {nSol=1; myAng = 0. ;}
01476 else if (myCos<=-1.) {nSol=1; myAng = M_PI ;}
01477 else {nSol=2; myAng = acos(myCos) ;}
01478 s = ( myAng -(alfa-beta))/Rho1;
01479 if (s<0) s+= 2.*M_PI/fabs(Rho1);
01480 s1 = s;
01481 if (nSol==2) {
01482 s =(-myAng -(alfa-beta))/Rho1;
01483 if (s< 0) s+= 2.*M_PI/fabs(Rho1);
01484 if (s<s1) s1 = s;
01485 }
01486 } else {
01487 TComplex A = CP1/CP2*(Im*Rho2);
01488 TComplex B = 1.-CdX/CP2*(Im*Rho2);
01489 double cba[3]={B.Rho2()-1., 2*(A.Re()*B.Re()+A.Im()*B.Im()), A.Rho2()};
01490 nSol = SqEqu(cba, L);
01491 if (nSol< 0) break;
01492 if (nSol==0) nSol=1;
01493 if (L[0]>0) s1=L[0];
01494 if (nSol>1 && L[1]>0 && L[1] <s1) s1 = L[1];
01495 }
01496
01497 break;
01498 }
01499 case 0: {
01500 if (fabs((CdX/CP1).Im())>fabs((CP2/CP1).Im())*1e6) break;
01501 nSol = 1;
01502 s = (CdX/CP2).Im()/(CP1/CP2).Im();
01503 if (s<0) break;
01504 s1 = s;
01505 break;
01506 }
01507 default: assert(0);
01508 }
01509 if (nSol) {
01510 TCircle tc1(*one),tc2(*two);
01511 double xy[2];
01512 tc1.Eval(s1,xy);
01513 s = tc2.Path(xy);
01514 if (s<0 && kase) s += 2.*M_PI/fabs(Rho2);
01515 if (s>0 ) s2 = s;
01516 }
01517
01518
01519 if (one != this) {s=s1;s1=s2;s2=s;}
01520 if (S2) *S2=s2;
01521 return s1;
01522 }
01523
01524 void TCircle::Test4()
01525 {
01526 double pars[4][2][5] = {
01527 {{0,0,1,0.,-0.001},{0,0.0,1,0,0.002}},
01528 {{0,0,1,0.,-0.001},{0,0.1,1,0,0.002}},
01529 {{0,0,1,0.,-0.001},{0,0.0,1,1,1e-8 }},
01530 {{0,0,1,0.,-1e-8 },{0,0.0,1,1,1e-8 }}};
01531 for (int ip=0;ip<4;ip++) {
01532 TCircle tc1(pars[ip][0],pars[ip][0]+2,pars[ip][0][4]);
01533 TCircle tc2(pars[ip][1],pars[ip][1]+2,pars[ip][1][4]);
01534 tc1.Move(-20);
01535 tc2.Move(-20);
01536 double s2;
01537 double s1 = tc1.Path(tc2,&s2);
01538 printf("s1=%g s2=%g \n",s1,s2);
01539 }
01540 }
01541
01542 double TCircle::Eval(double step,double *X, double *D) const
01543 {
01544
01545 sgCX1 =TComplex(fX[0],fX[1]);
01546 sgCD1 =TComplex(fD[0],fD[1]);
01547 sgImTet =TComplex(0,step*fRho);
01548 sgCOne =expOne(sgImTet);
01549 sgCf1 =sgImTet*sgCOne;
01550
01551 sgCD2 = sgCD1*sgCf1+sgCD1;
01552 sgCX2 = sgCD1*sgCOne*step;
01553 if (X) {
01554 X[0] = sgCX2.Re()+sgCX1.Re();
01555 X[1] = sgCX2.Im()+sgCX1.Im();
01556 }
01557 if (D) {
01558 sgCD2/= TComplex::Abs(sgCD2);
01559 D[0] = sgCD2.Re();
01560 D[1] = sgCD2.Im();
01561 }
01562 return step;
01563 }
01564
01565 double TCircle::Move(double step)
01566 {
01567 Eval(step,fX,fD);
01568 if (fEmx && fEmx->mHH>0) MoveErrs(step);
01569 if (fabs(fD[0])>1) fD[0]=(fD[0]<0)? -1:1;
01570 if (fabs(fD[1])>1) fD[1]=(fD[1]<0)? -1:1;
01571 return step;
01572 }
01573
01574 void TCircle::MakeMtx(double S,double F[3][3])
01575 {
01576 enum {kH=0,kA,kC};
01577 memset(F[0],0,sizeof(F[0][0])*3*3);
01578 F[kH][kH] = sgCf1.Re()+1.;
01579 double dSdH = sgCf1.Im();
01580
01581 F[kH][kA] = S*sgCOne.Re();
01582 double dSdA = S*sgCOne.Im();
01583 TComplex llCOneD = S*S*expOneD(-sgImTet);
01584 F[kH][kC] = llCOneD.Re();
01585 double dSdC = llCOneD.Im();
01586
01587 F[kA][kH] = -dSdH*fRho;
01588 F[kA][kA] = 1-dSdA*fRho;
01589 F[kA][kC] = S+dSdC*fRho;
01590 F[kC][kC] = 1;
01591 }
01592
01593 void TCircle::MoveErrs(double s)
01594 {
01595 double F[3][3];
01596 MakeMtx(s,F);
01597 fEmx->Move(F);
01598 }
01599
01600 void TCircle::Rot(double angle)
01601 {
01602 Rot(cos(angle),sin(angle));
01603 }
01604
01605 void TCircle::Rot(double cosa,double sina)
01606 {
01607 TComplex CX(fX[0],fX[1]),CP(fD[0],fD[1]);
01608 TComplex A (cosa,sina);
01609 CX *=A; fX[0] = CX.Re(); fX[1]=CX.Im();
01610 CP *=A; CP/=TComplex::Abs(CP);
01611 fD[0] = CP.Re(); fD[1]=CP.Im();
01612 }
01613
01614 void TCircle::Backward()
01615 {
01616 fRho=-fRho;fD[0]=-fD[0];fD[1]=-fD[1];
01617 if(fEmx) fEmx->Backward();
01618 }
01619
01620
01621 void TCircle::Test2()
01622 {
01623
01624
01625
01626
01627
01628
01629
01630
01631
01632
01633
01634
01635
01636
01637
01638
01639 }
01640
01641 void TCircle::Test3()
01642 {
01643
01644
01645
01646
01647
01648
01649
01650
01651
01652
01653
01654
01655
01656
01657
01658
01659
01660
01661
01662
01663
01664
01665
01666
01667
01668
01669
01670
01671
01672 }
01673
01674 void TCircle::TestMtx()
01675 {
01676 double Dir[8],X[8]={0},Rho[2],step,F[3][3],Del[3],Dif[3]={0};
01677 double maxEps = 0;
01678 int nErr=0,tally=0;
01679 for (int iR = 1010;iR>-1010;iR-=20) {
01680 Rho[0] = 1./iR;
01681 int iAlf = 360*gRandom->Rndm();
01682 Del[0] = 0.001;
01683 Del[1] = M_PI/180*0.01;
01684 Del[2] = 1e-4+ Rho[0]*0.001;
01685 for (int iStep=10;iStep<=350;iStep+=10){
01686 Dir[0] = cos(iAlf/180.*M_PI);
01687 Dir[1] = sin(iAlf/180.*M_PI);
01688 TCircle tc(X,Dir,Rho[0]);
01689 step = iStep/180.*M_PI*abs(iR);
01690 tc.Eval(step,X+2,Dir+2);
01691 tc.MakeMtx(step,F);
01692
01693 for (int iHAR=0;iHAR<3;iHAR++) {
01694 double minFak = 1e-4;
01695 for (double fak=1;fak>minFak;fak/=2) {
01696
01697 memcpy(X +4,X ,sizeof(X[0] )*2);
01698 memcpy(Dir+4,Dir,sizeof(Dir[0])*2);
01699 Rho[1]=Rho[0];
01700 switch (iHAR) {
01701 case 0: {
01702 X[4+0] = X[0]-Dir[1]*Del[0]*fak;
01703 X[4+1] = X[1]+Dir[0]*Del[0]*fak;
01704 break;}
01705
01706 case 1: {
01707 Dir[4+0] = cos(iAlf/180.*M_PI+Del[1]*fak);
01708 Dir[4+1] = sin(iAlf/180.*M_PI+Del[1]*fak);
01709 break;}
01710
01711 case 2: {
01712 Rho[1] = Rho[0]+Del[2]*fak;
01713 break;}
01714 }
01715
01716 TCircle tcc(X+4,Dir+4,Rho[1]);
01717 tcc.Move(step);
01718 double myStep = tcc.Path(X+2);
01719 tcc.Move(myStep);
01720 tcc.Eval(0,X+6,Dir+6);
01721 double dX[2];
01722 TCL::vsub(X+6,X+2,dX,2);
01723
01724 myStep = -TCL::vdot(dX,Dir+2,2)/TCL::vdot(Dir+6,Dir+2,2);
01725 tcc.Eval(myStep,X+6,Dir+6);
01726
01727 for (int jHAR=0;jHAR<2; jHAR++) {
01728 switch(jHAR) {
01729 case 0: {
01730 Dif[0] = (X[6+0]-X[2+0])*(-Dir[2+1])
01731 + (X[6+1]-X[2+1])*( Dir[2+0]);
01732 break;}
01733 case 1: {
01734 Dif[1] = (atan2(Dir[6+1],Dir[6+0])
01735 - atan2(Dir[2+1],Dir[2+0]));
01736 if (Dif[1]> M_PI) Dif[1]-=2*M_PI;
01737 if (Dif[1]< -M_PI) Dif[1]+=2*M_PI;
01738 break;}
01739 }
01740 tally++;
01741 double mat = F[jHAR][iHAR];
01742
01743
01744
01745
01746 double est = Dif[jHAR]/(Del[iHAR]*fak);
01747 double min = Del[jHAR]/Del[iHAR]*0.001;
01748
01749
01750 if (jHAR==iHAR) {
01751 if (fabs(est) < 1e-4) est = 0;
01752 if (fabs(mat) < 1e-4) mat = 0;
01753 }
01754
01755 double eps =(fabs(est-mat)*2)
01756 /(fabs(est)+fabs(mat)+min);
01757 if (eps < 1e-2) {
01758 if (eps>maxEps) maxEps=eps;
01759 break;
01760 }
01761 if (fak > minFak*2) continue;
01762 nErr++;
01763 if (eps>maxEps) maxEps=eps;
01764 printf("%6d Mtx[%d][%d] \t%g \t%g \tAngle=%d \tRad=%d \tLen=%g\n",
01765 tally,jHAR,iHAR,F[jHAR][iHAR],est,
01766 iAlf,iR,myStep);
01767
01768 } } } } }
01769 printf("TestMtx: %d errors maxEps=%g\n",nErr,maxEps);
01770
01771 }
01772
01773
01774
01775
01776 TCircleFitter::TCircleFitter()
01777 {
01778 Clear();
01779 SetEmx();
01780 }
01781
01782 void TCircleFitter::Clear(const char*)
01783 {
01784 fArr.Reset();
01785 memset(fBeg,0,fEnd-fBeg+1);
01786 TCircle::Clear();
01787 }
01788
01789 TCircleFitterAux* TCircleFitter::GetAux(int i) const
01790 {
01791 return (TCircleFitterAux*)(fArr.GetArray()+i*TCircleFitterAux::dSize());
01792 }
01793
01794 const double* TCircleFitter::GetX(int i) const
01795 {
01796 return &(fAux[i].x);
01797 }
01798
01799 double* TCircleFitter::GetX(int i)
01800 {
01801 return &(fAux[i].x);
01802 }
01803
01804 void TCircleFitter::Add(double x,double y,const double *errs)
01805 {
01806 fNuse =++fN;
01807 int n = fN*TCircleFitterAux::dSize();
01808 if (fArr.GetSize()<n) {fArr.Set(n*2);fAux=0;}
01809 if (!fAux) fAux = GetAux(0);
01810 TCircleFitterAux *aux = fAux+fN-1;
01811 aux->x = x; aux->y=y; aux->exy[0]=1; aux->exy[2]=1; aux->ezz=1;aux->wt=0;
01812 if (errs) AddErr(errs);
01813 }
01814
01815 void TCircleFitter::Add(double x,double y,double z)
01816 {
01817 fNuse =++fN;
01818 int n = fN*TCircleFitterAux::dSize();
01819 if (fArr.GetSize()<n) {fArr.Set(n*2);fAux=0;}
01820 if (!fAux) fAux = GetAux(0);
01821 TCircleFitterAux *aux = fAux+fN-1;
01822 aux->x = x; aux->y=y; aux->z=z;
01823 aux->exy[0]=1; aux->exy[1]=0; aux->exy[2]=1;aux->ezz=1;aux->wt=0;
01824 }
01825
01826 void TCircleFitter::AddErr(const double *errs,double ezz)
01827 {
01828 TCircleFitterAux *aux = fAux+fN-1;
01829 assert(ezz>=0);
01830 assert(errs[0]>=0);
01831 assert(errs[2]>=0);
01832 double spur = errs[0]+errs[2];
01833 assert(spur>0);
01834 double *e = aux->exy;
01835 memcpy(e,errs,sizeof(aux->exy));
01836
01837 if (e[0]<0 && e[0]>-1e-5*spur) {e[0]=0;e[1]=0;}
01838 if (e[2]<0 && e[2]>-1e-5*spur) {e[2]=0;e[1]=0;}
01839 assert(e[1]*e[1]<=1.01*e[0]*e[2]);
01840
01841 aux->wt = 0;
01842 aux->ezz = ezz;
01843 }
01844
01845 void TCircleFitter::AddErr(double errhh,double ezz)
01846 {
01847 TCircleFitterAux *aux = fAux+fN-1;
01848 assert(ezz>=0);
01849 assert(errhh>0);
01850 aux->wt = 1./errhh;
01851 aux->exy[0]= -1;
01852 aux->ezz = ezz;
01853 }
01854
01855 void TCircleFitter::AddZ(double z,double ez)
01856 {
01857
01858 fAux[fN-1].z =z;
01859 fAux[fN-1].ezz=ez;
01860 }
01861
01862 double TCircleFitter::Fit()
01863 {
01864 static const int nAVERs = &fRr-&fXx;
01865 static int nCall=0; nCall++;
01866 double xx, yy, xx2, yy2;
01867 double f, g, h, p, q, t, g0, g02, d=0;
01868 double xroot, ff, fp;
01869 double dx, dy, xnom,wt,hord,tmp,radius2,radiuc2;
01870 fKase = fCase;
01871 if (fNuse < 3) return 3e33;
01872 TCircleFitterAux *aux = GetAux(0);
01873 dx = aux[fN-1].x - aux[0].x;
01874 dy = aux[fN-1].y - aux[0].y;
01875 hord = sqrt(dx*dx+dy*dy);
01876 fCos = dx/hord;
01877 fSin = dy/hord;
01878 fNor[0] = -fSin,fNor[1] = fCos;
01879 const double *exy=0;
01880 for (int iter=0;iter<2;iter++) {
01881 fWtot = 0;
01882 memset(&fXgravity,0,sizeof(double)*(nAVERs+2));
01883 for (int i=0; i<fN; i++) {
01884 if (aux[i].wt<0) continue;
01885 int kase = iter;
01886 if (aux[i].wt >0) kase+=2;
01887 if (aux[i].exy[0]>0) kase+=4;
01888 switch (kase) {
01889 case 0: wt = 1; break;
01890 case 2: ; case 3: wt = aux[i].wt; break;
01891 case 5: ; case 7: {
01892 fNor[0] = fXCenter - aux[i].x;
01893 fNor[1] = fYCenter - aux[i].y;
01894 tmp = sqrt(fNor[0]*fNor[0]+fNor[1]*fNor[1]);
01895 fNor[0]/=tmp; fNor[1]/=tmp; }
01896 case 4:; case 6: {
01897 exy = aux[i].exy;
01898 wt = (fNor[0]*fNor[0]*exy[0]
01899 +fNor[0]*fNor[1]*exy[1]*2
01900 +fNor[1]*fNor[1]*exy[2]);
01901 if (wt<1e-8) wt = 1e-8;
01902 wt = 1/wt; break;}
01903 default: assert(0);
01904 }
01905 aux[i].wt = wt;
01906 fWtot += wt;
01907 fXgravity += aux[i].x *wt;
01908 fYgravity += aux[i].y *wt;
01909 }
01910 fXgravity /= fWtot;
01911 fYgravity /= fWtot;
01912
01913 for (int i=0; i<fN; i++) {
01914 dx = aux[i].x-fXgravity;
01915 dy = aux[i].y-fYgravity;
01916 xx = dx*fCos + dy*fSin;
01917 yy = -dx*fSin + dy*fCos;
01918 xx2 = xx*xx;
01919 yy2 = yy*yy;
01920 wt = aux[i].wt;
01921 fXx += xx2 *wt;
01922 fYy += yy2 *wt;
01923 fXy += xx*yy *wt;
01924 fXrr += xx*(xx2+yy2) *wt;
01925 fYrr += yy*(xx2+yy2) *wt;
01926 fRrrr += (xx2+yy2)*(xx2+yy2) *wt;
01927 }
01928 double *dd = &fXx;
01929 for (int i=0;i<nAVERs;i++) {dd[i]/=fWtot;}
01930 fRr = fXx+fYy;
01931
01932 if (fNuse <= 3 && !fKase) fKase=1;
01933 if (!fKase) fKase =(fYy < fXx *(0.5)*(0.5)/210)? 1:2;
01934 SWIT: switch(fKase) {
01935 case 1: {
01936
01937
01938
01939
01940
01941
01942
01943 double myCof[3];
01944 fPol[0] = 1;
01945 fPol[1] = 0; fPol[2] = 1./(2*sqrt(fXx));
01946 fPol[3] = fRr; fPol[4] = fXrr/(2*fXx); fPol[5] = 1.;
01947 double tmp = sqrt(fPol[3]*fPol[3]
01948 +fPol[4]*fPol[4]*(4*fXx )
01949 +fPol[5]*fPol[5]*(fRrrr )
01950 +fPol[3]*fPol[5]*(-fRr ) *2
01951 +fPol[4]*fPol[5]*(-2*fXrr) *2);
01952 fPol[3]/=tmp;fPol[4]/=tmp;fPol[5]/=tmp;
01953 myCof[0] = 0;
01954 myCof[1] = - (fPol[2]*(4*fXy));
01955 myCof[2] = - (fPol[4]*(4*fXy) + fPol[5]*(-2*fYrr));
01956 fC = myCof[0]*fPol[0]+myCof[1]*fPol[1]+myCof[2]*fPol[3];
01957 fA = myCof[1]*fPol[2]+myCof[2]*fPol[4];
01958 fB = myCof[2]*fPol[5];
01959 fYd = (fabs(fB)>1e-6) ? 1./fB : 1e6;
01960 fXd = fA*fYd;
01961 }
01962 break;
01963
01964 case 2: {
01965
01966 f = (3.*fXx+fYy);
01967 g = (fXx+3.*fYy);
01968 h = 2*fXy;
01969 p = fXrr;
01970 q = fYrr;
01971 t = fRrrr;
01972 g0 = (fXx+fYy);
01973 g02 = g0*g0;
01974 fA = -4.0;
01975 fB = (f*g-t-h*h)/g02;
01976 fC = (t*(f+g)-2.*(p*p+q*q))/(g02*g0);
01977 d = (t*(h*h-f*g)+2.*(p*p*g+q*q*f)-4.*p*q*h)/(g02*g02);
01978 xroot = 1.0;
01979 for (int i=0; i<5; i++) {
01980 ff = (((xroot+fA)*xroot+fB)*xroot+fC)*xroot+d;
01981 fp = ((4.*xroot+3.*fA)*xroot+2.*fB)*xroot+fC;
01982 xroot -= ff/fp;
01983 }
01984 fG1 = xroot*g0;
01985 xnom = (g-fG1)*(f-fG1)-h*h;
01986
01987
01988 if (xnom<1e-20) { fKase=1; goto SWIT;}
01989
01990 fXd = ( p*(g-fG1)-q*h )/xnom;
01991 fYd = (-p*h +q*(f-fG1))/xnom;
01992 }
01993 break;
01994
01995 default: assert(0);
01996 }
01997 fXCenter = fXd*fCos-fYd*fSin + fXgravity;
01998 fYCenter = fXd*fSin+fYd*fCos + fYgravity;
01999 if (!exy) break;
02000 }
02001
02002
02003 switch (fKase) {
02004 case 1: {
02005 fCorrR = sqrt(1+fA*fA+fC*fB );
02006 fCorrB = sqrt(1+fA*fA );
02007 fRho = fabs(fB)/fCorrR;
02008 int sgB = (fB<0)? -1:1;
02009 fNy = sgB/sqrt(1+fA*fA);
02010 fNx = fA*fNy;
02011 fH = -fC*sgB/(fCorrR+fCorrB);
02012 fChi2 = (4*fA*fXy +4*fYy- 2*fB*fYrr)/4;
02013 fChi2 /= (fCorrR*fCorrR);
02014 fX[0] = fNx*fH; fX[1] = fNy*fH;
02015
02016
02017 fD[0] = fNy; fD[1] =-fNx;
02018
02019 Rot(fCos,fSin);
02020 fX[0] += fXgravity;
02021 fX[1] += fYgravity;
02022 }
02023 break;
02024 case 2: {
02025 radiuc2 = fXd*fXd+fYd*fYd;
02026 radius2 = radiuc2+fG1;
02027 fR = ::sqrt(radius2);
02028 fRd = ::sqrt(radiuc2);
02029 fRho = 1./fR;
02030 fH = -fG1/(fR+fRd);
02031 if (fRd>1e-6) { fNx = fXd/fRd; fNy = fYd/fRd;}
02032 else { fNx = 0; fNy = 1; fRd = 0; }
02033 fChi2 = (fG1-fRr)/2;
02034 fX[0] = fNx*fH; fX[1] = fNy*fH;
02035
02036
02037 fD[0] = fNy; fD[1] =-fNx;
02038
02039 Rot(fCos,fSin);
02040 fX[0] += fXgravity;
02041 fX[1] += fYgravity;
02042
02043
02044
02045
02046
02047 }
02048 break;
02049 default: assert(0);
02050 }
02051 fNdf = fNuse-3;
02052 if (fNdf>0) fChi2 *= fWtot/fNdf;
02053 tmp = fD[0]*(aux[0].x-fX[0])+fD[1]*(aux[0].y-fX[1]);
02054
02055 fBack = 0;
02056 if (tmp>0) {fD[0]*=-1;fD[1]*=-1;fRho*=-1;fBack=1;}
02057 return fChi2;
02058 }
02059
02060 void TCircleFitter::MakeErrs()
02061 {
02062 fEmx->Clear();
02063 double F[3][3]; memset(F[0],0,sizeof(F));
02064 double myFact = 1.;
02065 switch (fKase) {
02066 case 1: {
02067 fCov[0] = fPol[2]*fPol[2]+ fPol[4]*fPol[4];
02068 fCov[1] = fPol[4]*fPol[5];
02069 fCov[2] = fPol[5]*fPol[5];
02070 fCov[3] = fPol[1]*fPol[2]+ fPol[3]*fPol[4];
02071 fCov[4] = fPol[3]*fPol[5];
02072 fCov[5] = fPol[0]*fPol[0]+ fPol[1]*fPol[1]+fPol[3]*fPol[3];
02073 for (int i=0;i<6;i++) {fCov[i]*=4;}
02074 int sgB = (fB<0)? -1:1;
02075 double corrRB = fCorrR+fCorrB;
02076 double corrR3 = fCorrR*fCorrR*fCorrR;
02077
02078
02079 F[0][0] = sgB*fA*fC/(corrRB*fCorrB*fCorrR);
02080 F[0][1] = 0.5*sgB*fC*fC/(corrRB*corrRB*fCorrR);
02081 F[0][2] = 0.5*sgB*fC*fB/(corrRB*corrRB*fCorrR)
02082 - sgB /(corrRB );
02083 F[1][0] = -1/(fCorrB*fCorrB);
02084 F[2][0] = - sgB*fA*fB/(corrR3);
02085 F[2][1] = -0.5*sgB*fC*fB/(corrR3)+sgB/fCorrR;
02086 F[2][2] = -0.5*sgB*fB*fB/(corrR3);
02087 myFact = (fCorrR*fCorrR);
02088 break;
02089 }
02090 case 2: {
02091
02092 double aRho = fabs(fRho),aRho2 = aRho*aRho,aRho3 = aRho*aRho2;
02093 for (int i=0,li=0;i< 3;li+=++i) {
02094 for (int j=0 ;j<=i;j++ ) {
02095 fCov[li+j] = d2F(i,j)*0.5; }}
02096 TCL::trsinv(fCov,fCov ,3);
02097 #if 1
02098
02099 double dRho = -fH/(fRd*fR);
02100 dRho = 1/fRd-1/fR;
02101 F[0][0] = fXd*dRho;
02102 F[0][1] = fYd*dRho;
02103 F[0][2] = -0.5*aRho;
02104 F[1][0] = -fNy*aRho;
02105 F[1][1] = fNx*aRho;
02106 F[1][2] = 0;
02107 F[2][0] = -fXd*aRho3;
02108 F[2][1] = -fYd*aRho3;
02109 F[2][2] = -0.5*aRho3;
02110 #endif
02111 #if 0
02112 fXCenter = fXd*fCos-fYd*fSin + fXgravity;
02113 fYCenter = fXd*fSin+fYd*fCos + fYgravity;
02114 fNx = fXd*fCos-fYd*fSin;
02115 fNy = fXd*fSin+fYd*fCos;
02116 double nor = sqrt(fNx*fNx+fNy*fNy);
02117 fNx/=nor; fNy/=nor;
02118
02119
02120
02121 fX[0] = fXCenter -fNx*fR;
02122 fX[1] = fYCenter -fNy*fR;
02123 fD[0] = fNy;
02124 fD[1] =-fNx;
02125
02126
02127 F[0][0] = ( fNx*fCos+fNy*fSin) -fXd*aRho;
02128 F[0][1] = (-fNx*fSin+fNy*fCos) -fYd*aRho;
02129 F[0][2] = -0.5*aRho;
02130 F[1][0] = (-fNy*fCos + fNx*fSin)*aRho;
02131 F[1][1] = ( fNy*fSin + fNx*fCos)*aRho;
02132 F[1][2] = 0;
02133 F[2][0] = -fXd*aRho3;
02134 F[2][1] = -fYd*aRho3;
02135 F[2][2] = -0.5*aRho3;
02136
02137 #endif
02138
02139 break;
02140 }
02141 default: assert(0);
02142 }
02143 TCL::vscale(fCov,myFact/fWtot,fCov,6);
02144 TCL::trasat(F[0],fCov,fEmx->Arr(),3,3);
02145 if (fBack) {fEmx->mHA*=-1;fEmx->mAC*=-1;}
02146 }
02147
02148 double TCircleFitter::EvalChi2()
02149 {
02150 if (!fNuse) return 0;
02151 TCircle M(this);
02152 double sum = 0,wtot=0,wt;
02153 TCircleFitterAux *aux = GetAux(0);
02154 const double *p = M.Pos();
02155 for (int i=0;i<fN;i++) {
02156 if (aux[i].wt<0) continue;
02157 double s = M.Path(&(aux[i].x));
02158 M.Move(s);
02159 wt = aux[i].wt;
02160 sum+= (pow(p[0]-aux[i].x,2)+pow(p[1]-aux[i].y,2))*wt;
02161 wtot +=wt;
02162 }
02163 if (fNdf) sum /= fNdf;
02164 fChi2 = sum;
02165 return sum;
02166 }
02167
02168 double TCircleFitter::FixAt(const double vals[5],int flag)
02169 {
02177
02178 assert(fEmx);
02179 assert(flag);
02180 double g[6]={1,0,1,0,0,1},c[6]={1,0,1,0,0,1}
02181 ,e[6],adj[3]={0},amda[3],dlt[2];
02182 int sel[3] ={!!(flag&1), !!(flag&2), !!(flag&4)};
02183 int nFix=0;
02184 if (sel[0]) {
02185 nFix++;
02186 dlt[0] = vals[0]-fX[0]; dlt[1] = vals[1]-fX[1];
02187 adj[0] = -dlt[0]*fD[1]+dlt[1]*fD[0];
02188 }
02189 if (sel[1]) {
02190 nFix++;
02191 adj[1] = vals[3]-atan2(fD[1],fD[0]);
02192 if (adj[1]< -M_PI) adj[1] += 2*M_PI;
02193 if (adj[1]> M_PI) adj[1] -= 2*M_PI;
02194 }
02195 if (sel[2]) {
02196 nFix++;
02197 adj[2] = vals[4]-fRho;
02198 }
02199
02200 for (int i=0,li=0;i< 3;li+=++i) {
02201 for (int j=0 ;j<=i;j++ ) {
02202 if (!(sel[i]&sel[j])) continue;
02203 c[li+j] = (*fEmx)[li+j]; } }
02204 double addXi2=0;
02205 TCL::trsinv(c ,c ,3 );
02206 TCL::trasat(adj,c,&addXi2,1,3);
02207
02208 TCL::trsinv(fEmx->Arr(),e,3);
02209 for (int i=0,li=0;i< 3;li+=++i) {
02210 for (int j=0 ;j<=i;j++ ) {
02211 if (!(sel[i]|sel[j])) continue;
02212 e[li+j] = (i==j);
02213 if (!(sel[i]&sel[j])) continue;
02214 g[li+j] = (*fEmx)[li+j];
02215 } }
02216 TCL::trsinv(g ,g ,3 );
02217 TCL::trsa (g ,adj ,amda,3,1);
02218 TCL::trsa (fEmx->Arr(),amda,adj ,3,1);
02219 TCL::trsinv(e ,fEmx->Arr(),3 );
02220
02221 for (int i=0,li=0;i< 3;li+=++i) {if (sel[i]) (*fEmx)[li+i]=0;}
02222
02223 fX[0] += -adj[0]*fD[1];
02224 fX[1] += adj[0]*fD[0];
02225
02226
02227
02228
02229 double S = sin(adj[1]);
02230 double C = cos(adj[1]);
02231 double d0 = fD[0];
02232 fD[0] = d0*C-fD[1]*S;
02233 fD[1] = d0*S+fD[1]*C;
02234
02235 fRho += adj[2];
02236 fNdf+=nFix;
02237 fChi2 += (addXi2-fChi2*nFix)/fNdf;
02238 return fChi2;
02239 }
02240
02241 void TCircleFitter::Skip(int idx)
02242 {
02243 fAux[idx].exy[0] = -1;
02244 SetNdf(fNdf-1);
02245 }
02246
02247 void TCircleFitter::SetNdf(int ndf)
02248 {
02249 fChi2*=fNdf; if (ndf) fChi2/=ndf; fNdf=ndf;
02250 }
02251
02252 void TCircleFitter::Print(const char* txt) const
02253 {
02254 if (!txt) txt="";
02255 printf("TCircleFitter::NPoints = %d method=%d",fN,fKase);
02256 if (fChi2) printf(" Chi2 = %g",fChi2);
02257 printf("\n");
02258 TCircle::Print();
02259
02260 int iP = (strchr(txt,'P') || strchr(txt,'p'));
02261 int iE = (strchr(txt,'E') || strchr(txt,'e'));
02262 int iF = (strchr(txt,'F') || strchr(txt,'f'));
02263 int iZ = (strchr(txt,'Z') || strchr(txt,'z'));if(iZ){};
02264 TCircleFitterAux *aux = GetAux(0);
02265 if (iP) {
02266 for (int ip=0;ip<fN;ip++) {
02267 printf("%3d - X: %g\tY:%g \tZ:%g",ip,aux[ip].x,aux[ip].y,aux[ip].z);
02268 if (iE)
02269 printf(" \tExy: %g %g %g \tEz:%g "
02270 ,aux[ip].exy[0],aux[ip].exy[1],aux[ip].exy[2],aux[ip].ezz);
02271 printf("\n");
02272 }}
02273 if (iF) {
02274 TCircle circ(this);
02275 const double *xy = GetX(0);
02276 double ds=circ.Path(xy);
02277 circ.Move(ds);
02278 double s=0;
02279 for (int i=0;i<fN;i++) {
02280 xy = GetX(i);
02281 ds = circ.Path(xy);
02282 s+=ds;
02283 circ.Move(ds);
02284 if (fabs( s)<1e-6) s=0;
02285 if (fabs(ds)<1e-6)ds=0;
02286 printf("%3d - S=%g(%g) \tX: %g=%g \tY:%g=%g \tdirX=%g dirY=%g\n"
02287 ,i,s,ds
02288 ,xy[0],circ.Pos()[0]
02289 ,xy[1],circ.Pos()[1]
02290 ,circ.Dir()[0],circ.Dir()[1]);
02291 }}
02292
02293 }
02294
02295 void TCircleFitter::Test(int iTest)
02296 {
02297
02298 int myKase = (iTest/1 )%10;
02299 int negRad = (iTest/10 )%10;
02300 int bakPts = (iTest/100 )%10;
02301 int noShft = (iTest/1000)%10;
02302
02303
02304 enum {nPts=5};
02305 double e[4],x[3],dir[2];
02306 double aShift[6];
02307 aShift[0]=-acos(0.25);
02308 aShift[1]=-acos(0.50);
02309 aShift[2]= 0;
02310 aShift[3]= acos(0.25);
02311 aShift[5]= acos(0.50);
02312 memset(aShift,0,sizeof(aShift));
02313 double xCenter[2];
02314 double RERR = 0.01;
02315 TRandom ran;
02316
02317 static TCanvas* myCanvas[9]={0};
02318 static TH1F *hh[10]={0};
02319 static const char *hNams[]={"pH","pA","pC","Xi2","dErr"};
02320 static const double lims[][2]={{-5 ,5},{-5 ,5 },{-5 ,5}
02321 ,{ 0 ,5},{-0.05,0.05}};
02322 const int nPads = sizeof(hNams)/sizeof(void*);
02323
02324 if(!myCanvas[0]) myCanvas[0]=new TCanvas("TCircleFitter_Test0","",600,800);
02325 myCanvas[0]->Clear();
02326 myCanvas[0]->Divide(1,nPads);
02327
02328 for (int i=0;i<nPads;i++) {
02329 delete hh[i]; hh[i]= new TH1F(hNams[i],hNams[i],100,lims[i][0],lims[i][1]);
02330 myCanvas[0]->cd(i+1); hh[i]->Draw();
02331 }
02332
02333 static TH1F *hc[10]={0};
02334 static const char *cNams[]={"HHpull","HApull","AApull","HCpull","ACpull","CCpull"};
02335 static const int cPads=sizeof(cNams)/sizeof(char*);
02336 if(!myCanvas[1]) myCanvas[1]=new TCanvas("TCircleFitter_TestCorr1","",600,800);
02337 myCanvas[1]->Clear();
02338 myCanvas[1]->Divide(1,cPads);
02339
02340 for (int i=0;i<cPads;i++) {
02341 delete hc[i]; hc[i]= new TH1F(cNams[i],cNams[i],100,-15,15);
02342 myCanvas[1]->cd(i+1); hc[i]->Draw();
02343 }
02344
02345 static TH1F *hl[10]={0};
02346 static const char *lNams[]={"XP","YP","GP","XY","XG","YG"};
02347 static const int lPads=sizeof(lNams)/sizeof(char*);
02348 if(!myCanvas[2]) myCanvas[2]=new TCanvas("TCircleFitter_TestCorr2","",600,800);
02349 myCanvas[2]->Clear();
02350 myCanvas[2]->Divide(1,lPads);
02351
02352 for (int i=0;i<lPads;i++) {
02353 delete hl[i]; hl[i]= new TH1F(lNams[i],lNams[i],100,-5,5);
02354 myCanvas[2]->cd(i+1); hl[i]->Draw();
02355 }
02356
02357
02358 int nFit = 0,isgn;
02359 static uint seed=0;
02360
02361
02362 for (int ir = 50; ir <= 50; ir +=5) {
02363 double aR = ir;
02364 double len = 0.3*aR*3;
02365 for (double ang0 = 0; ang0 <= 0; ang0+=0.05) {
02366 double R = (negRad)? -aR:aR;
02367 double dang = len/R/(nPts-1);
02368 if (bakPts) dang*=-1;
02369 double C0 = cos(ang0);
02370 double S0 = sin(ang0);
02371
02372 double myX[2]={0,0},myD[2]={C0,S0};
02373 xCenter[0] = myX[0]-myD[1]*R;
02374 xCenter[1] = myX[1]+myD[0]*R;
02375 double corr[6]={0},korr[6]={0};
02376 seed++;
02377 for (int iter=0;iter<2;iter++) {
02378 static const int nEv = 100000;
02379 ran.SetSeed(seed);
02380 for (int ev=0;ev<nEv;ev++) {
02381 TCircleFitter helx;
02382 TCircle idex(myX,myD,1/R);;
02383 for (int is=0;is<nPts;is++) {
02384 double ang = ang0 + dang*is;
02385 double S = sin(ang),C = cos(ang);
02386 double eR = ran.Gaus(0,RERR);
02387 double shift = aShift[is%5];
02388 shift=0;
02389 double SS = sin(ang+shift);
02390 double CC = cos(ang+shift);
02391 e[0] = pow(RERR*SS,2);
02392 e[1] =-pow(RERR ,2)*CC*SS;
02393 e[2] = pow(RERR*CC,2);
02394
02395 x[0] = myX[0] + (R)*(S-S0);
02396 x[1] = myX[1] - (R)*(C-C0);
02397 x[0]+= -SS*eR;
02398 x[1]+= CC*eR;
02399
02400 helx.Add (x[0],x[1]);
02401 helx.AddErr (RERR*RERR);
02402 }
02403 helx.SetCase(myKase);
02404 double Xi2 = helx.Fit();
02405 double Xj2 = helx.EvalChi2();
02406 assert(fabs(Xi2-Xj2) < 1e-2*(Xi2+Xj2+0.01));
02407
02408 assert (!myKase || helx.GetCase()==myKase);
02409
02410 nFit++;
02411 helx.MakeErrs();
02412 if ((isgn=helx.Emx()->Sign())<0) {
02413 ::Warning("Test1","Negative errmtx %d",isgn);continue;}
02414
02415 if (helx.Rho()*idex.Rho() < 0) idex.Backward();
02416
02417 double myShift = (aR*M_PI/180)*360*gRandom->Rndm();
02418 if (noShft) myShift=0;
02419 helx.Move(myShift);
02420 idex.Move(myShift);
02421 if ((isgn=helx.Emx()->Sign())<0) {
02422 ::Warning("Test2","Negative errmtx %d",isgn);continue;}
02423 double s = idex.Path(helx.Pos());
02424 idex.Eval(s,x,dir);
02425
02426 double dd[10];
02427 double dx = helx.Pos()[0]-x[0];
02428 double dy = helx.Pos()[1]-x[1];
02429 dd[0] = -dx*dir[1]+dy*dir[0];
02430 dd[1] = atan2(helx.Dir()[1],helx.Dir()[0])-atan2(dir[1],dir[0]);
02431 if (dd[1]> M_PI) dd[1]-=2*M_PI;
02432 if (dd[1]<-M_PI) dd[1]+=2*M_PI;
02433 dd[2] = helx.Rho()-idex.Rho();
02434
02435 if (!iter) {
02436 for (int i=0,li=0;i< 3;li+=++i) {
02437 for (int j=0;j<=i;j++) {
02438 corr[li+j]+= (*(helx.Emx()))[li+j];
02439 korr[li+j]+= dd[i]*dd[j];
02440 } }
02441 continue;
02442 }
02443
02444 dd[3] = dd[0]/sqrt(corr[0]);
02445 dd[4] = dd[1]/sqrt(corr[2]);
02446 dd[5] = dd[2]/sqrt(corr[5]);
02447 dd[6] = Xi2;
02448 dd[7] = sqrt(helx.Emx()->mHH)-RERR;
02449
02450 for (int ih=0;ih<nPads;ih++) { hh[ih]->Fill(dd[ih+3]);}
02451
02452
02453 double cc[6]={0},dia[3];
02454 for (int i=0,li=0;i< 3;li+=++i) {
02455 dia[i] = corr[li+i];
02456 for (int j=0;j<=i;j++) {
02457 cc[li+j] = (dd[i]*dd[j]-corr[li+j])/sqrt(dia[i]*dia[j]);
02458 }
02459 }
02460
02461 for (int ih=0;ih<cPads;ih++) { hc[ih]->Fill(cc[ih]);}
02462
02463 double myCenter[3];
02464 myCenter[0]=xCenter[0]-helx.fXgravity;
02465 myCenter[1]=xCenter[1]-helx.fYgravity;
02466 double xx = myCenter[0];
02467 myCenter[0] = xx*helx.fCos + myCenter[1]*helx.fSin;
02468 myCenter[1] = -xx*helx.fSin + myCenter[1]*helx.fCos;
02469 myCenter[2] = R*R - (pow(myCenter[0],2)+pow(myCenter[1],2));
02470 for (int j=0;j<3;j++) {dd[j]=(&(helx.fXd))[j]-myCenter[j];}
02471 dd[3+0] = dd[0]/sqrt(helx.fCov[0]);
02472 dd[3+1] = dd[1]/sqrt(helx.fCov[2]);
02473 dd[3+2] = dd[2]/sqrt(helx.fCov[5]);
02474 for (int ih=0;ih<3;ih++) { hl[ih]->Fill(dd[ih+3]);}
02475 cc[0] = (dd[3+0]*dd[3+1]-helx.fCov[1])/sqrt(helx.fCov[0]*helx.fCov[2]);
02476 cc[1] = (dd[3+0]*dd[3+2]-helx.fCov[3])/sqrt(helx.fCov[0]*helx.fCov[5]);
02477 cc[2] = (dd[3+1]*dd[3+2]-helx.fCov[4])/sqrt(helx.fCov[2]*helx.fCov[5]);
02478 for (int ih=0;ih<3;ih++) { hl[ih+3]->Fill(cc[ih]);}
02479 }
02480 if (iter) break;
02481 TCL::vscale(corr,1./nEv,corr,6);
02482 TCL::vscale(korr,1./nEv,korr,6);
02483 double dia[3];
02484 for (int i=0,li=0;i< 3;li+=++i) {
02485 dia[i]=sqrt(corr[li+i]);
02486 int n = 0;
02487 for (int j=0;j<=i;j++) {n+=printf("%15g",corr[li+j]/(dia[i]*dia[j]));}
02488 n = 45-n;
02489 for (int j=0;j< n;j++) {printf(" ");}
02490 for (int j=0;j<=i;j++) {n+=printf("%15g",korr[li+j]/(dia[i]*dia[j]));}
02491 printf("\n");
02492 }
02493
02494
02495 }
02496 }
02497 }
02498
02499 for (int i=0;myCanvas[i];i++) {
02500 myCanvas[i]->Modified();myCanvas[i]->Update();}
02501
02502 while(!gSystem->ProcessEvents()){gSystem->Sleep(200);};
02503
02504 }
02505
02506
02507 void TCircleFitter::TestCorr(int kase)
02508 {
02509
02510
02511
02512
02513
02514
02515 if (!(kase&3 ))kase+=1+2;
02516 if (!(kase&12))kase+=4+8;
02517 enum {nPts=20};
02518 double e[4],x[3],ex[3];
02519 double aShift[6];
02520 aShift[0]=-acos(0.25);
02521 aShift[1]=-acos(0.50);
02522 aShift[2]= 0;
02523 aShift[3]= acos(0.25);
02524 aShift[5]= acos(0.50);
02525 double RERR = 0.001;
02526 TRandom ran;
02527 static TCanvas* myCanvas=0;
02528 static TH1F *hh[6]={0,0,0,0,0,0};
02529 static const char *hNams[]={"HA","HA-","HC","HC-","AC","AC-",0};
02530 if(!myCanvas) myCanvas=new TCanvas("TCircleFitter_TestCorr","",600,800);
02531 myCanvas->Clear();
02532 myCanvas->Divide(1,6);
02533
02534 for (int i=0;i<6;i++) {
02535 delete hh[i]; hh[i]= new TH1F(hNams[i],hNams[i],100,-1,1);
02536 myCanvas->cd(i+1); hh[i]->Draw();
02537 }
02538
02539 int nFit = 0;
02540 for (int ir = 50; ir <= 1000; ir +=5) {
02541 double aR = ir;
02542 double len = 100; if (len>aR*3) len = aR*3;
02543 for (double ang0 = -3; ang0 < 3.1; ang0+=0.05) {
02544 for (int sgn = -1; sgn<=1; sgn+=2) {
02545 if ((sgn>0) && !(kase&4)) continue;
02546 if ((sgn<0) && !(kase&8)) continue;
02547 double R = sgn*aR;
02548 double dang = len/R/nPts;
02549 double C0 = cos(ang0);
02550 double S0 = sin(ang0);
02551 TCircleFitter helx;
02552 for (int is=0;is<nPts;is++) {
02553 double ang = ang0 + dang*is;
02554 double S = sin(ang),C = cos(ang);
02555 double eR = ran.Gaus(0,RERR)*sgn;
02556 double shift = aShift[is%5];
02557
02558 double SS = sin(ang+shift);
02559 double CC = cos(ang+shift);
02560 e[0] = pow(RERR*SS,2);
02561 e[1] =-pow(RERR ,2)*CC*SS;
02562 e[2] = pow(RERR*CC,2);
02563
02564 x[0] = 100 + (R)*(S-S0);
02565 x[1] = 200 - (R)*(C-C0);
02566 ex[0]= x[0]-SS*eR;
02567 ex[1]= x[1]+CC*eR;
02568 helx.Add (ex[0],ex[1],e);
02569 }
02570 helx.Fit();
02571 if (!(helx.GetCase()&kase)) continue;
02572 nFit++;
02573 helx.MakeErrs();
02574 int iFix = 0;
02575 if (kase&16) iFix +=1;
02576 if (kase&32) iFix +=4;
02577 if (iFix) {
02578 double vals[5];
02579 TCL::ucopy(x,vals,3);
02580 vals[3]=0;
02581 vals[4]=1./R;
02582 helx.FixAt(vals,iFix);
02583 }
02584 x[0] = 100 ;
02585 x[1] = 200 ;
02586 double s = helx.Path(x);
02587 assert(s<0);
02588 assert(fabs(s) < len);
02589 helx.Move(s);
02590 double dd[6],hf[6];
02591 double dx = helx.Pos()[0]-x[0];
02592 double dy = helx.Pos()[1]-x[1];
02593 const TCEmx_t *emx = helx.Emx();
02594 dd[0] = -dx*S0+dy*C0;
02595 dd[1] = atan2(helx.Dir()[1],helx.Dir()[0])-ang0;
02596 if (dd[1]> M_PI) dd[1]-=2*M_PI;
02597 if (dd[1]<-M_PI) dd[1]+=2*M_PI;
02598 dd[2] = helx.Rho()-1./R;
02599 hf[0] = (dd[0]*dd[1]) *1e1/(RERR*RERR);
02600 hf[1] = (emx->mHA) *1e1/(RERR*RERR);
02601 hf[2] = dd[0]*dd[2] *1e3/(RERR*RERR);
02602 hf[3] = (emx->mHC) *1e3/(RERR*RERR);
02603 hf[4] = dd[1]*dd[2] *1e4/(RERR*RERR);
02604 hf[5] = (emx->mAC) *1e4/(RERR*RERR);
02605
02606
02607 for (int ih=0;ih<6;ih++) { hh[ih]->Fill(hf[ih]);}
02608 }
02609 }
02610 }
02611 myCanvas->Modified();
02612 myCanvas->Update();
02613 while(!gSystem->ProcessEvents()){gSystem->Sleep(200);};
02614
02615 }
02616
02617 void TCircle::Show(int nPts,const double *Pts,int pstep)
02618 {
02619 static TCanvas *myCanvas = 0;
02620 static TGraph *ptGraph = 0;
02621 static TGraph *ciGraph = 0;
02622 double x[100],y[100];
02623 if (nPts>100) nPts=100;
02624 for (int i=0;i<nPts;i++) { x[i]=Pts[i*pstep+0]; y[i]=Pts[i*pstep+1]; }
02625
02626
02627 if(!myCanvas) myCanvas = new TCanvas("TCircle_Show","",600,800);
02628 myCanvas->Clear();
02629 delete ptGraph; delete ciGraph;
02630
02631 ptGraph = new TGraph(nPts , x, y);
02632 ptGraph->SetMarkerColor(kRed);
02633 ptGraph->Draw("A*");
02634
02635 TCircle tc(this);
02636 double xy[2];
02637 xy[0] = x[0];
02638 xy[1] = y[0];
02639 double s = tc.Path(xy);
02640 tc.Move(s);
02641 xy[0] = x[nPts-1];
02642 xy[1] = y[nPts-1];
02643 s = tc.Path(xy);
02644 if (s<0) { tc.Backward(); s = tc.Path(xy);}
02645 double ds = s/99;
02646 for (int i=0;i<100;i++) {x[i]=tc.Pos()[0];y[i]=tc.Pos()[1];tc.Move(ds);}
02647
02648 ciGraph = new TGraph(100 , x, y);
02649 ciGraph->Draw("Same CP");
02650 myCanvas->Modified();
02651 myCanvas->Update();
02652 while(!gSystem->ProcessEvents()){gSystem->Sleep(200);};
02653
02654 }
02655
02656 THelixFitter::THelixFitter():fPoli1Fitter(1)
02657 {
02658 Clear();
02659 SetEmx();
02660 }
02661
02662 THelixFitter::~THelixFitter()
02663 {;}
02664
02665 void THelixFitter::Clear(const char*)
02666 {
02667 fCircleFitter.Clear();
02668 fPoli1Fitter.Clear();
02669 fPoli1Fitter.SetCoefs(1);
02670 fChi2=0;
02671 }
02672
02673 void THelixFitter::Print(const char*) const
02674 {
02675 THelixTrack::Print();
02676 fCircleFitter.Print();
02677 fPoli1Fitter.Print();
02678 }
02679
02680 void THelixFitter::Add (double x,double y,double z)
02681 {
02682 fCircleFitter.Add(x,y,z);
02683 }
02684
02685 void THelixFitter::AddErr(const double *err2xy,double err2z)
02686 {
02687 fCircleFitter.AddErr(err2xy,err2z);
02688 }
02689
02690 void THelixFitter::AddErr(double errhh,double err2z)
02691 {
02692 fCircleFitter.AddErr(errhh,err2z);
02693 }
02694
02695 double THelixFitter::Fit()
02696 {
02697 TCircleFitterAux* myAux= GetAux(0);
02698 int nDat = Size();
02699 double Xi2xy = fCircleFitter.Fit();
02700 if (Xi2xy>1e11) return Xi2xy;
02701 int ndfXY = fCircleFitter.Ndf();
02702 TCircle circ(fCircleFitter);
02703 const double *xy=0;
02704 xy = &(myAux[nDat-1].x);
02705 double z1 = xy[2];
02706 double s1 = circ.Path(xy);
02707
02708 xy = &(myAux[0].x);
02709 double z0 = xy[2];
02710 double s = circ.Path(xy);
02711
02712 double tanDip = (z1-z0)/(s1-s);
02713
02714 circ.Move(s);
02715
02716 for (int iDat=0;iDat<nDat;iDat++) {
02717 TCircleFitterAux* aux = myAux+iDat;
02718 xy = &(aux->x);
02719 double ds = circ.Path(xy,aux->exy);
02720 circ.Move(ds); s+=ds;
02721
02722 double corErr = 0;;
02723 if (aux->exy[0]>0) {
02724 const double *dc = circ.Dir();
02725 corErr = tanDip*tanDip*
02726 (dc[0]*dc[0]*aux->exy[0]
02727 +dc[1]*dc[1]*aux->exy[2]
02728 +dc[0]*dc[1]*aux->exy[1]*2);
02729 }
02730 fPoli1Fitter.Add(s,aux->z,aux->ezz+corErr);
02731 }
02732 double Xi2z = fPoli1Fitter.Fit();
02733
02734 int ndfSz = fPoli1Fitter.Ndf();
02735 Update(1);
02736 int ndf = ndfSz+ndfXY;
02737 fChi2 = Xi2xy*ndfXY+Xi2z*ndfSz;
02738 if (ndf) fChi2/=ndf;
02739 return fChi2;
02740
02741 }
02742
02743 double THelixFitter::FixAt(const double val[5],int flag)
02744 {
02745 double xx[3],s;
02746 memcpy(xx,fX,sizeof(xx));
02747 int move = (flag&1);
02748 if (move) {
02749 s = fCircleFitter.Path(val);
02750 fCircleFitter.Move(s);
02751 fPoli1Fitter.Move(s);
02752 }
02753 double Xi2c = fCircleFitter.FixAt(val,flag);
02754 if (flag&1) fPoli1Fitter.FixAt(0.,val[2]);
02755
02756 if (move) {
02757 s = fCircleFitter.Path(xx);
02758 fCircleFitter.Move(s);
02759 fPoli1Fitter.Move(s);
02760 }
02761 Update(1+2);
02762
02763 double Xi2z = fPoli1Fitter.Chi2();
02764 int ndfc = fCircleFitter.Ndf();
02765 int ndfz = fPoli1Fitter.Ndf();
02766
02767 int ndf = ndfc+ndfz;
02768 fChi2 = Xi2c*ndfc+Xi2z*ndfz;
02769 if (ndf) fChi2/=ndf;
02770 return fChi2;
02771 }
02772
02773 void THelixFitter::Skip(int idx)
02774 {
02775 fCircleFitter.Skip(idx);
02776 fPoli1Fitter.Skip(idx);
02777 int ndfc = fCircleFitter.Ndf();
02778 int ndfz = fPoli1Fitter.Ndf();
02779 int ndf = ndfc+ndfz;
02780 fChi2 = fCircleFitter.Chi2()*ndfc+fPoli1Fitter.Chi2()*ndfz;
02781 if (ndf) fChi2/=ndf;
02782 }
02783
02784 void THelixFitter::Update(int kase)
02785 {
02786 if(kase&1) {
02787 const double *pol = fPoli1Fitter.Coe();
02788 fCosL = 1./sqrt(pol[1]*pol[1]+1);
02789 double *haslet = fCircleFitter.Pos();
02790 fX[0] = haslet[0];
02791 fX[1] = haslet[1];
02792 fX[2] = pol[0];
02793 fP[0] = haslet[2]*fCosL;
02794 fP[1] = haslet[3]*fCosL;
02795 fP[2] = pol[1]*fCosL;
02796 fRho = haslet[4];
02797 }
02798 if(kase&2) {
02799 double emx[3];
02800 emx[0] = fPoli1Fitter.Emx()[0];
02801 emx[1] = fPoli1Fitter.Emx()[1]*fCosL*fCosL;
02802 emx[2] = fPoli1Fitter.Emx()[2]*fCosL*fCosL*fCosL*fCosL;
02803 fEmx->Set(fCircleFitter.Emx()->Arr(),emx);
02804 }
02805 }
02806
02807 void THelixFitter::MakeErrs()
02808 {
02809 fCircleFitter.MakeErrs();
02810 fPoli1Fitter.MakeErrs();
02811 Update(2);
02812 }
02813
02814 double THelixFitter::EvalChi2()
02815 {
02816 double Xi2c = fCircleFitter.EvalChi2();
02817 double Xi2z = fPoli1Fitter.EvalChi2();
02818 fChi2 = Xi2c*fCircleFitter.Ndf()+Xi2z*fPoli1Fitter.Ndf();
02819 fChi2/=(fCircleFitter.Ndf()+fPoli1Fitter.Ndf()+1e-10);
02820 return fChi2;
02821 }
02822
02823 void THelixFitter::Test(int kase)
02824 {
02825
02826
02827
02828
02829
02830
02831
02832
02833 if (!(kase&3 ))kase+=1+2;
02834 if (!(kase&12))kase+=4+8;
02835
02836 enum {nPts=5,nHH=8};
02837 double e[4],x[3],xe[3];
02838 double aShift[6];
02839 aShift[0]=-acos(0.25);
02840 aShift[1]=-acos(0.50);
02841 aShift[2]= 0;
02842 aShift[3]= acos(0.25);
02843 aShift[5]= acos(0.50);
02844 double RERR = 0.1;
02845 double ZERR = 0.1;
02846 TRandom ran;
02847 static TCanvas* myCanvas[9]={0};
02848 static TH1F *hh[nHH]={0};
02849 static const char *hNams[]={"pH","pA","pC","pZ","pD","Xi2","Xi2E","Xi2d",0};
02850 if(!myCanvas[0]) myCanvas[0]=new TCanvas("THelixFitter_TestC1","",600,800);
02851 myCanvas[0]->Clear();
02852 myCanvas[0]->Divide(1,nHH);
02853
02854 for (int i=0;i<nHH;i++) {
02855 double low = (i>=5)? 0:-5;
02856 double upp = 5;
02857 delete hh[i]; hh[i]= new TH1F(hNams[i],hNams[i],100,low,upp);
02858 myCanvas[0]->cd(i+1); hh[i]->Draw();
02859 }
02860
02861
02862 static TH1F *h2h[4]={0,0,0,0};
02863 static const char *h2Nams[]={"targYY","targZZ","targYZ","calcYZ",0};
02864 int n2h=4;
02865 if(!myCanvas[1]) myCanvas[1]=new TCanvas("THelixFitter_TestC2","",600,800);
02866 myCanvas[1]->Clear();
02867 myCanvas[1]->Divide(1,n2h);
02868 for (int i=0;i<n2h;i++) {
02869 delete h2h[i]; h2h[i]= new TH1F(h2Nams[i],h2Nams[i],100,-5,5);
02870 myCanvas[1]->cd(i+1); h2h[i]->Draw();
02871 }
02872
02873
02874
02875 static TH1F *h3h[4]={0,0,0,0};
02876 static const char *h3Nams[]={"dcaXY","dcaXYNor","dcaZ","dcaZNor",0};
02877 int n3h=4;
02878 if(!myCanvas[2]) myCanvas[2]=new TCanvas("THelixFitter_TestC3","",600,800);
02879 myCanvas[2]->Clear();
02880 myCanvas[2]->Divide(1,n3h);
02881 for (int i=0;i<n3h;i++) {
02882 delete h3h[i]; h3h[i]= new TH1F(h3Nams[i],h3Nams[i],100,-5,5);
02883 myCanvas[2]->cd(i+1); h3h[i]->Draw();
02884 }
02885
02886
02887
02888 double spotSurf[4]= {-100,1,0,0};
02889 double spotAxis[3][3]= {{0,1,0},{0,0,1},{1,0,0}};
02890
02891
02892 int nFit = 0,isgn;
02893 for (double idip=-1;idip<=1;idip+=0.2){
02894 double dip = idip;
02895
02896 double cosDip = cos(dip);
02897 double sinDip = sin(dip);
02898 double tanDip = tan(dip); if(tanDip){};
02899 for (int ir = 50; ir <= 1000; ir +=20) {
02900 double aR = ir;
02901 double len = 100; if (len>aR*3) len = aR*3;
02902 for (double ang00 = -3; ang00 < 3.1; ang00+=0.2) {
02903 double ang0 = ang00;
02904
02905 for (int sgn = -1; sgn<=1; sgn+=2) {
02906 if(sgn>0 && !(kase&4)) continue;
02907 if(sgn<0 && !(kase&8)) continue;
02908
02909 double R = sgn*aR;
02910 double dang = len/R/nPts;
02911 double C0 = cos(ang0);
02912 double S0 = sin(ang0);
02913 THelixFitter helx;
02914
02915 double trakPars[7]={100,200,300,C0*cosDip,S0*cosDip,sinDip,1/R};
02916 THelixTrack trak(trakPars+0,trakPars+3,trakPars[6]);
02917
02918 for (int is=0;is<nPts;is++) {
02919 double ang = ang0 + dang*is;
02920 double S = sin(ang),C = cos(ang);
02921 double eR = ran.Gaus(0,RERR)*sgn;
02922 double eZ = ran.Gaus(0,ZERR);
02923 double shift = aShift[is%5];
02924
02925 double SS = sin(ang+shift);
02926 double CC = cos(ang+shift);
02927 e[0] = pow(RERR*SS,2);
02928 e[1] =-pow(RERR ,2)*CC*SS;
02929 e[2] = pow(RERR*CC,2);
02930 e[3] = pow(ZERR,2);
02931 x[0] = 100 + (R)*(S-S0);
02932 x[1] = 200 - (R)*(C-C0);
02933 double len = (R)*(ang-ang0);
02934 x[2] = 300 + len*tan(dip);
02935 xe[0]= x[0]-SS*eR;
02936 xe[1]= x[1]+CC*eR;
02937 xe[2]= x[2]+eZ;
02938 helx.Add (xe[0],xe[1],xe[2]);
02939 helx.AddErr(e,e[3]);
02940 }
02941 double Xi2 =helx.Fit();
02942 if(!(kase&helx.GetCase())) continue;
02943
02944 helx.MakeErrs();
02945 if ((isgn=helx.Emx()->Sign())<0) {
02946 ::Warning("Test1","Negative errmtx %d",isgn);continue;}
02947 nFit++;
02948 if (kase&16) Xi2=helx.FixAt(x);
02949
02950
02951 if (kase&32) {
02952 double vals[5];
02953 TCL::ucopy(x,vals,3);vals[3]=0;vals[4]=1./R;
02954 Xi2=helx.FixAt(vals,4);
02955 }
02956 if (kase&128) helx.Show();
02957 double Xi2E =helx.EvalChi2();
02958
02959 trak.Move(0.3*len/cosDip);
02960 memcpy(x,trak.Pos(),sizeof(x));
02961 ang0 = atan2(trak.Dir()[1],trak.Dir()[0]);
02962
02963 double s = helx.Path(x);
02964
02965
02966
02967 double pos[3],dir[3],rho;
02968 helx.Move(s);
02969 if ((isgn=helx.Emx()->Sign())<0) {
02970 ::Warning("Test2","Negative errmtx %d",isgn);continue;}
02971 THEmx_t *emx = helx.Emx();
02972 helx.Get (pos,dir,rho);
02973 double psi = atan2(dir[1],dir[0]);
02974 double sinPsi = sin(psi);
02975 double cosPsi = cos(psi);
02976 double tanPsi = sinPsi/cosPsi; if(tanPsi){};
02977 double dd[10],hf[10];
02978 double dx = x[0]-pos[0];
02979 double dy = x[1]-pos[1];
02980 dd[0] = -dx*sinPsi+dy*cosPsi;
02981 hf[0] = dd[0]/sqrt(emx->mHH+1e-20);
02982 dd[2] = psi-ang0;
02983 if (dd[2]> M_PI) dd[2]-=2*M_PI;
02984 if (dd[2]<-M_PI) dd[2]+=2*M_PI;
02985 hf[1] = dd[2]/sqrt(emx->mAA+1e-20);
02986 dd[4] = rho-1./R;
02987 hf[2] = dd[4]/sqrt(emx->mCC+1e-20);
02988 dd[6] = (helx.Pos()[2]-x[2])/pow(helx.GetCos(),2);
02989 hf[3] = dd[6]/sqrt(emx->mZZ+1e-20);
02990 dd[8] = asin(dir[2])-dip;
02991 if (dd[8]> M_PI) dd[8]-=2*M_PI;
02992 if (dd[8]<-M_PI) dd[8]+=2*M_PI;
02993 hf[4] = dd[8]/(sqrt(emx->mLL));
02994 hf[5] = Xi2;
02995 hf[6] = Xi2E;
02996 hf[7] = Xi2E-Xi2+1;
02997 for (int ih=0;ih<nHH;ih++) { hh[ih]->Fill(hf[ih]);}
02998
02999
03000 double xIde[3],pIde[3],xFit[3],pFit[3],eSpot[3],hfil,sIde,sFit;
03001
03002 int closePoint=0;
03003 spotSurf[0] = -110;
03004
03005 { spotSurf[0] = -x[0]; closePoint=2006;}
03006 sIde = trak.Step(200.,spotSurf,4, xIde,pIde,closePoint);
03007
03008 if (fabs(spotSurf[0]+TCL::vdot(xIde,spotSurf+1,3))>0.001) {
03009 printf("***Wrong point found**\n");
03010 trak.Print();
03011 assert(0);
03012 }
03013 sFit = helx.Step(200.,spotSurf,4, xFit,pFit,closePoint);
03014 if (sFit>=1000 ) continue;
03015 if (fabs(pIde[0]-pFit[0])>0.1) continue;
03016 helx.Move(sFit);
03017 emx = helx.Emx();
03018 helx.GetSpot(spotAxis,eSpot);
03019 hfil = (xFit[1]-xIde[1]); hfil/= sqrt(eSpot[0]);
03020 h2h[0]->Fill(hfil);
03021 hfil = (xFit[2]-xIde[2]); hfil/= sqrt(eSpot[2]);
03022 h2h[1]->Fill(hfil);
03023 hfil = (xFit[1]-xIde[1])*(xFit[2]-xIde[2]);
03024 h2h[2]->Fill(hfil*100);
03025 h2h[3]->Fill(hfil/sqrt(eSpot[0]*eSpot[2]));
03026
03027
03028
03029
03030 double dcaXY,dcaZ,dcaEmx[3];
03031 double sDca = helx.Dca(trakPars,dcaXY,dcaZ,dcaEmx);
03032 if (fabs(sDca)<1000) {
03033 h3h[0]->Fill(dcaXY);
03034 h3h[1]->Fill(dcaXY/sqrt(dcaEmx[0]));
03035 h3h[2]->Fill(dcaZ );
03036 h3h[3]->Fill(dcaZ /sqrt(dcaEmx[2]));
03037 }
03038
03039
03040 }
03041 }
03042 }
03043 }
03044 for (int ih=0;myCanvas[ih];ih++) {
03045 myCanvas[ih]->Modified();
03046 myCanvas[ih]->Update();
03047 }
03048 while(!gSystem->ProcessEvents()){gSystem->Sleep(200);};
03049 }
03050
03051 void THelixFitter::Show() const
03052 {
03053 static TCanvas *myCanvas = 0;
03054 static TGraph *ptGraph[2] = {0,0};
03055 static TGraph *ciGraph[2] = {0,0};
03056 double x[100],y[100],z[100],l[100]
03057 , X[100],Y[100],Z[100];
03058 int nPts = Size();
03059 if (nPts>100) nPts=100;
03060 TCircleFitterAux* aux=GetAux(0);
03061 THelixTrack tc(this);
03062 double s = tc.Path(aux[0].x,aux[0].y); tc.Move(s);
03063 s = tc.Path(aux[nPts-1].x,aux[nPts-1].y);
03064 if (s<0) { tc.Backward();}
03065 l[0]=0;
03066 double ds=0;
03067 for (int i=0;i<nPts;i++) {
03068 if (i) {ds = tc.Path(aux[i].x,aux[i].y);tc.Move(ds);l[i]=l[i-1]+ds;}
03069 x[i]=aux[i].x; y[i]=aux[i].y; z[i]=aux[i].z;
03070 X[i]=tc.Pos()[0];Y[i]=tc.Pos()[1];Z[i]=tc.Pos()[2];
03071 }
03072
03073
03074 if(!myCanvas) myCanvas = new TCanvas("THelixFitter_Show","",600,800);
03075 myCanvas->Clear();
03076 myCanvas->Divide(1,2);
03077
03078 delete ptGraph[0]; delete ciGraph[0];
03079 ptGraph[0] = new TGraph(nPts , x, y);
03080 ptGraph[0]->SetMarkerColor(kRed);
03081 myCanvas->cd(1); ptGraph[0]->Draw("A*");
03082 delete ptGraph[1]; delete ciGraph[1];
03083 ptGraph[1] = new TGraph(nPts , l, z);
03084 ptGraph[1]->SetMarkerColor(kRed);
03085 myCanvas->cd(2); ptGraph[1]->Draw("A*");
03086
03087 ciGraph[0] = new TGraph(nPts , X, Y);
03088 myCanvas->cd(1); ciGraph[0]->Draw("Same CP");
03089 ciGraph[1] = new TGraph(nPts , l, Z);
03090 myCanvas->cd(2); ciGraph[1]->Draw("Same CP");
03091
03092 myCanvas->Modified();
03093 myCanvas->Update();
03094 while(!gSystem->ProcessEvents()){gSystem->Sleep(200);};
03095
03096 }
03097
03098 double TCircleFitter::f()
03099 {
03100 return 4*((fG1-fRr)/2)*fR*fR;
03101 }
03102
03103 double TCircleFitter::F()
03104 {
03105 return (fG1-fRr)/2;
03106 }
03107
03108 double TCircleFitter::df(int i)
03109 {
03110 switch (i) {
03111 case 0: return -4*(fXrr -2*fXx*fXd - 2*fXy*fYd);
03112 case 1: return -4*(fYrr -2*fXy*fXd - 2*fYy*fYd);
03113 case 2: return 2*(fG1 - fRr);
03114 default: assert(0);
03115 }
03116 assert(0);
03117 return 0;
03118 }
03119
03120 double TCircleFitter::d2f(int i,int j)
03121 {
03122
03123
03124
03125 assert(j<=i);
03126 int ij = i+10*j;
03127 switch(ij) {
03128 case 0: return 8*fXx;
03129 case 01: return 8*fXy;
03130 case 11: return 8*fYy;
03131 case 02:; case 12: return 0;
03132 case 22: return 2;
03133 default: printf ("Kase=%d\n",ij); assert(0);
03134 assert(0);
03135 }
03136 return 0;
03137 }
03138
03139 double TCircleFitter::Rho2 () { return fRho*fRho;}
03140
03141
03142 double TCircleFitter::dRho2(int i)
03143 {
03144 double ans =fRho*fRho;
03145 ans *= -ans;
03146 switch (i) {
03147 case 0: return 2*ans*fXd;
03148 case 1: return 2*ans*fYd;
03149 case 2: return ans;
03150 }
03151 assert(0);
03152 return 0;
03153 }
03154
03155 double TCircleFitter::d2Rho2(int i,int j)
03156 {
03157
03158
03159
03160 if (j>i) { int jj = j; j = i; i = jj;}
03161 int ij = i+10*j;
03162 double rho2 = fRho*fRho,rho4 = rho2*rho2, rho6 = rho4*rho2;
03163 switch(ij) {
03164 case 0: return 8*(rho6)*fXd*fXd-2*(rho4);;
03165 case 01: return 8*(rho6)*fXd*fYd;
03166 case 11: return 8*(rho6)*fYd*fYd-2*(rho4);
03167 case 02: return 4*(rho6)*fXd;
03168 case 12: return 4*(rho6)*fYd;
03169 case 22: return 2*(rho6);
03170 default: printf ("Kase=%d\n",ij); assert(0);
03171 }
03172 }
03173
03174 double TCircleFitter::dF(int i)
03175 {
03176
03177
03178 double ans = 1./4*(df(i)*Rho2() + f()*dRho2(i));
03179 return ans;
03180 }
03181
03182 double TCircleFitter::d2F(int i,int j)
03183 {
03184
03185
03186 double ans = 1./4*(d2f(i,j)*Rho2() +df(j)*dRho2(i)
03187 +df (i) *dRho2(j)+f() *d2Rho2(i,j));
03188 return ans;
03189 }
03190
03191
03192
03193
03194
03195
03196
03197
03198
03199
03200
03201
03202
03203
03204
03205
03206
03207
03208
03209
03210
03211
03212
03213
03214
03215
03216
03217
03218
03219
03220
03221
03222
03223
03224
03225
03226
03227
03228
03229
03230
03231
03232
03233
03234
03235
03236
03237
03238
03239
03240
03241
03242
03243
03244
03245
03246
03247
03248
03249
03250
03251
03252
03253
03254
03255
03256
03257
03258
03259
03260
03261
03262
03263
03264
03265
03266
03267
03268
03269
03270
03271
03272
03273
03274
03275
03276
03277
03278
03279
03280
03281
03282
03283
03284
03285
03286
03287
03288
03289
03290
03291
03292
03293
03294
03295
03296
03297
03298
03299
03300
03301
03302
03303
03304
03305
03306
03307
03308
03309
03310
03311
03312
03313
03314
03315
03316
03317
03318
03319
03320
03321