00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016 #include "BDSQuadStepper.hh"
00017 #include "G4ThreeVector.hh"
00018 #include "G4TransportationManager.hh"
00019
00020 using std::max;
00021 extern G4double BDSLocalRadiusOfCurvature;
00022 extern G4int event_number;
00023
00024 const int DEBUG = 0;
00025
00026 BDSQuadStepper::BDSQuadStepper(G4Mag_EqRhs *EqRhs)
00027 : G4MagIntegratorStepper(EqRhs,6)
00028
00029 {
00030 fPtrMagEqOfMot = EqRhs;
00031 }
00032
00033
00034 void BDSQuadStepper::AdvanceHelix( const G4double yIn[],
00035 G4ThreeVector Bfld,
00036 G4double h,
00037 G4double yQuad[])
00038 {
00039 const G4double *pIn = yIn+3;
00040 G4ThreeVector GlobalR = G4ThreeVector( yIn[0], yIn[1], yIn[2]);
00041 G4ThreeVector GlobalP = G4ThreeVector( pIn[0], pIn[1], pIn[2]);
00042 G4ThreeVector InitMomDir = GlobalP.unit();
00043 G4double InitPMag = GlobalP.mag();
00044 G4double charge = (fPtrMagEqOfMot->FCof())/c_light;
00045 G4double kappa = - fPtrMagEqOfMot->FCof()*itsBGrad/InitPMag;
00046
00047 if (DEBUG) G4cout << "BDSQuadStepper: step= " << h/m << " m" << G4endl
00048 << " x= " << yIn[0]/m << "m" << G4endl
00049 << " y= " << yIn[1]/m << "m" << G4endl
00050 << " z= " << yIn[2]/m << "m" << G4endl
00051 << " px= " << yIn[3]/GeV << "GeV/c" << G4endl
00052 << " py= " << yIn[4]/GeV << "GeV/c" << G4endl
00053 << " pz= " << yIn[5]/GeV << "GeV/c" << G4endl
00054 << " q= " << charge/eplus << "e" << G4endl
00055 << " dBy/dx= " << itsBGrad/(tesla/m) << "T/m" << G4endl
00056 << " k= " << kappa/(1./(m*m)) << "m^-2" << G4endl
00057 << G4endl;
00058
00059 G4Navigator* QuadNavigator=
00060 G4TransportationManager::GetTransportationManager()->
00061 GetNavigatorForTracking();
00062
00063 G4double h2=h*h;
00064
00065 G4ThreeVector LocalR,LocalRp ;
00066
00067
00068
00069 if(fabs(kappa)<1.e-12)
00070 {
00071 G4ThreeVector positionMove = h * InitMomDir;
00072
00073 yQuad[0] = yIn[0] + positionMove.x();
00074 yQuad[1] = yIn[1] + positionMove.y();
00075 yQuad[2] = yIn[2] + positionMove.z();
00076
00077 yQuad[3] = GlobalP.x();
00078 yQuad[4] = GlobalP.y();
00079 yQuad[5] = GlobalP.z();
00080
00081 itsDist=0;
00082 }
00083 else
00084 {
00085 h2=h*h;
00086
00087 G4AffineTransform GlobalAffine=QuadNavigator->
00088 GetGlobalToLocalTransform();
00089
00090 G4ThreeVector LocalR=GlobalAffine.TransformPoint(GlobalR);
00091 G4ThreeVector LocalRp=GlobalAffine.TransformAxis(InitMomDir);
00092
00093 if (DEBUG) G4cout << "BDSQuadStepper: initial point in local coordinates:" << G4endl
00094 << " x= " << LocalR[0]/m << "m" << G4endl
00095 << " y= " << LocalR[1]/m << "m" << G4endl
00096 << " z= " << LocalR[2]/m << "m" << G4endl
00097 << " x'= " << LocalRp[0] << G4endl
00098 << " y'= " << LocalRp[1] << G4endl
00099 << " z'= " << LocalRp[2] << G4endl
00100 << G4endl;
00101
00102 G4double x0,xp,y0,yp,z0,zp;
00103 G4double x1,x1p,y1,y1p,z1,z1p;
00104 x0=LocalR.x();
00105 y0=LocalR.y();
00106 z0=LocalR.z();
00107 xp=LocalRp.x();
00108 yp=LocalRp.y();
00109 zp=LocalRp.z();
00110
00111
00112 G4ThreeVector LocalRpp;
00113 LocalRpp.setX(-zp*x0);
00114 LocalRpp.setY( zp*y0);
00115 LocalRpp.setZ( x0*xp - y0*yp);
00116
00117 LocalRpp*=kappa;
00118
00119
00120 G4double R_1 = LocalRpp.mag();
00121 if (DEBUG) G4cout << " curvature= " << R_1*m << "m^-1" << G4endl;
00122 if(R_1>0.)
00123 {
00124 G4double R=1./R_1;
00125
00126
00127 BDSLocalRadiusOfCurvature=R;
00128
00129
00130 itsDist= h2/(8*R);
00131
00132
00133 if((fabs(zp)>0.99)&&(fabs(kappa)<1.e-6))
00134 {
00135 G4double rootK=sqrt(fabs(kappa*zp));
00136 G4double rootKh=rootK*h*zp;
00137 G4double X11,X12,X21,X22;
00138 G4double Y11,Y12,Y21,Y22;
00139
00140 if (kappa>0)
00141 {
00142 X11= cos(rootKh);
00143 X12= sin(rootKh)/rootK;
00144 X21=-fabs(kappa)*X12;
00145 X22= X11;
00146
00147 Y11= cosh(rootKh);
00148 Y12= sinh(rootKh)/rootK;
00149 Y21= fabs(kappa)*Y12;
00150 Y22= Y11;
00151 }
00152 else
00153 {
00154 X11= cosh(rootKh);
00155 X12= sinh(rootKh)/rootK;
00156 X21= fabs(kappa)*X12;
00157 X22= X11;
00158
00159 Y11= cos(rootKh);
00160 Y12= sin(rootKh)/rootK;
00161 Y21= -fabs(kappa)*Y12;
00162 Y22= Y11;
00163 }
00164
00165 x1 = X11*x0 + X12*xp;
00166 x1p= X21*x0 + X22*xp;
00167
00168 y1 = Y11*y0 + Y12*yp;
00169 y1p= Y21*y0 + Y22*yp;
00170
00171 z1p=sqrt(1 - x1p*x1p -y1p*y1p);
00172
00173 G4double dx=x1-x0;
00174 G4double dy=y1-y0;
00175
00176
00177 G4double dR2=dx*dx+dy*dy;
00178 G4double dz=sqrt(h2*(1.-h2/(12*R*R))-dR2);
00179
00180
00181 G4double ScaleFac=(dx*dx+dy*dy+dz*dz)/h2;
00182 if(ScaleFac>1.0000001)
00183 {
00184 ScaleFac=sqrt(ScaleFac);
00185 dx/=ScaleFac;
00186 dy/=ScaleFac;
00187 dz/=ScaleFac;
00188 x1=x0+dx;
00189 y1=y0+dy;
00190 }
00191 z1=z0+dz;
00192 }
00193 else
00194
00195 {
00196
00197 G4double quadX= - kappa*x0*zp;
00198 G4double quadY= kappa*y0*zp;
00199 G4double quadZ= kappa*(x0*xp - y0*yp);
00200
00201
00202 G4double maxCurv=max(fabs(quadX),fabs(quadY));
00203 maxCurv=max(maxCurv,fabs(quadZ));
00204
00205 x1 = x0 + h*xp + quadX*h2/2;
00206 y1 = y0 + h*yp + quadY*h2/2;
00207 z1 = z0 + h*zp + quadZ*h2/2;
00208
00209 x1p = xp + quadX*h;
00210 y1p = yp + quadY*h;
00211 z1p = zp + quadZ*h;
00212
00213
00214 G4double quadX_end= - kappa*x1*z1p;
00215 G4double quadY_end= kappa*y1*z1p;
00216 G4double quadZ_end= kappa*(x1*x1p - y1*y1p);
00217
00218
00219 maxCurv=max(maxCurv,fabs(quadX_end));
00220 maxCurv=max(maxCurv,fabs(quadY_end));
00221 maxCurv=max(maxCurv,fabs(quadZ_end));
00222
00223 itsDist=maxCurv*h2/4.;
00224
00225
00226 G4double quadX_av=(quadX+quadX_end)/2;
00227 G4double quadY_av=(quadY+quadY_end)/2;
00228 G4double quadZ_av=(quadZ+quadZ_end)/2;
00229
00230 G4double x_prime_av=(xp + x1p)/2;
00231 G4double y_prime_av=(yp + y1p)/2;
00232 G4double z_prime_av=(zp + z1p)/2;
00233
00234 x1 = x0 + h*x_prime_av + quadX_av * h2/2;
00235 y1 = y0 + h*y_prime_av + quadY_av * h2/2;
00236 z1 = z0 + h*z_prime_av + quadZ_av * h2/2;
00237
00238 x1p = xp + quadX_av*h;
00239 y1p = yp + quadY_av*h;
00240 z1p = zp + quadZ_av*h;
00241
00242 G4double dx = (x1-x0);
00243 G4double dy = (y1-y0);
00244 G4double dz = (z1-z0);
00245 G4double chord2 = dx*dx + dy*dy + dz*dz;
00246 if(chord2>h2)
00247 {
00248 G4double hnew=h*sqrt(h2/chord2);
00249 x1=x0 + hnew*x_prime_av + quadX_av * hnew*hnew/2;
00250 y1=y0 + hnew*y_prime_av + quadY_av * hnew*hnew/2;
00251 z1=z0 + hnew*z_prime_av + quadZ_av * hnew*hnew/2;
00252
00253 x1p=xp + quadX_av*hnew;
00254 y1p=yp + quadY_av*hnew;
00255 z1p=zp + quadZ_av*hnew;
00256 }
00257 }
00258
00259 LocalR.setX(x1);
00260 LocalR.setY(y1);
00261 LocalR.setZ(z1);
00262
00263 LocalRp.setX(x1p);
00264 LocalRp.setY(y1p);
00265 LocalRp.setZ(z1p);
00266
00267 }
00268 else
00269 {
00270 LocalR += h*LocalRp;
00271 itsDist=0.;
00272 }
00273
00274 if (DEBUG) G4cout << "BDSQuadStepper: final point in local coordinates:" << G4endl
00275 << " x= " << LocalR[0]/m << "m" << G4endl
00276 << " y= " << LocalR[1]/m << "m" << G4endl
00277 << " z= " << LocalR[2]/m << "m" << G4endl
00278 << " x'= " << LocalRp[0] << G4endl
00279 << " y'= " << LocalRp[1] << G4endl
00280 << " z'= " << LocalRp[2] << G4endl
00281 << G4endl;
00282
00283
00284
00285 G4AffineTransform LocalAffine=QuadNavigator->
00286 GetLocalToGlobalTransform();
00287
00288 GlobalR = LocalAffine.TransformPoint(LocalR);
00289 GlobalP = InitPMag*LocalAffine.TransformAxis(LocalRp);
00290
00291 yQuad[0] = GlobalR.x();
00292 yQuad[1] = GlobalR.y();
00293 yQuad[2] = GlobalR.z();
00294
00295 yQuad[3] = GlobalP.x();
00296 yQuad[4] = GlobalP.y();
00297 yQuad[5] = GlobalP.z();
00298
00299 }
00300 }
00301
00302
00303 void BDSQuadStepper::Stepper( const G4double yInput[],
00304 const G4double dydx[],
00305 const G4double hstep,
00306 G4double yOut[],
00307 G4double yErr[] )
00308 {
00309 const G4int nvar = 6 ;
00310 G4int i;
00311
00312 const G4double *pIn = yInput+3;
00313 G4ThreeVector GlobalP = G4ThreeVector( pIn[0], pIn[1], pIn[2]);
00314 G4double InitPMag = GlobalP.mag();
00315 G4double kappa= - fPtrMagEqOfMot->FCof()*itsBGrad/InitPMag;
00316
00317 if(fabs(kappa)<1.e-6)
00318 {
00319 for(i=0;i<nvar;i++) yErr[i]=0;
00320 AdvanceHelix(yInput,0,hstep,yOut);
00321 }
00322 else
00323 {
00324 G4double yTemp[7], yIn[7];
00325
00326
00327
00328 for(i=0;i<nvar;i++) yIn[i]=yInput[i];
00329
00330
00331 G4double h = hstep * 0.5;
00332 AdvanceHelix(yIn, 0, h, yTemp);
00333 AdvanceHelix(yTemp, 0, h, yOut);
00334
00335
00336 h = hstep ;
00337 AdvanceHelix(yIn, 0, h, yTemp);
00338
00339 for(i=0;i<nvar;i++) yErr[i] = yOut[i] - yTemp[i] ;
00340 }
00341 return;
00342 }
00343
00344 G4double BDSQuadStepper::DistChord() const
00345 {
00346 return itsDist;
00347
00348
00349 }
00350
00351 BDSQuadStepper::~BDSQuadStepper()
00352 {}