/home/cern/BDSIM_new/src/myQuadStepper.cc

00001 /* BDSIM code.    Version 1.0
00002    Author: Grahame A. Blair, Royal Holloway, Univ. of London.
00003    Copyright (c) 2002 by G.A.Blair.  ALL RIGHTS RESERVED. 
00004 */
00005 
00006 #include "BDSGlobalConstants.hh" 
00007 #include "myQuadStepper.hh"
00008 #include "G4ThreeVector.hh"
00009 #include "G4LineSection.hh"
00010 #include "G4TransportationManager.hh"
00011 
00012 using std::max;
00013 extern G4double BDSLocalRadiusOfCurvature;
00014 extern G4int event_number;
00015 
00016 myQuadStepper::myQuadStepper(G4Mag_EqRhs *EqRhs)
00017   : G4MagIntegratorStepper(EqRhs,6)  // integrate over 6 variables only !!
00018                                        // position & velocity
00019 {
00020   fPtrMagEqOfMot = EqRhs;
00021 }
00022 
00023 
00024 void myQuadStepper::AdvanceHelix( const G4double  yIn[],
00025                                   G4ThreeVector Bfld,
00026                                   G4double  h,
00027                                   G4double  yOut[])
00028 {
00029   //G4cout<<"myQuadStepper: advancing through "<<h/m<<"  m "<<
00030   //"x="<<yIn[0] /m<<" y="<<yIn[1]/m<<" z="<<yIn[2]/m<<G4endl; 
00031 
00032   const G4double *pIn = yIn+3;
00033   G4ThreeVector v0= G4ThreeVector( pIn[0], pIn[1], pIn[2]);  
00034 
00035   G4ThreeVector GlobalPosition= G4ThreeVector( yIn[0], yIn[1], yIn[2]);  
00036   G4double InitMag=v0.mag();
00037   G4ThreeVector InitMomDir=v0.unit();
00038   
00039   //G4double h2=h*h;
00040 
00041   G4Navigator* HelixNavigator=
00042     G4TransportationManager::GetTransportationManager()->
00043     GetNavigatorForTracking();
00044   
00045   G4AffineTransform LocalAffine=HelixNavigator-> 
00046     GetLocalToGlobalTransform();
00047   
00048   G4AffineTransform GlobalAffine=HelixNavigator->
00049     GetGlobalToLocalTransform();
00050 
00051   G4ThreeVector LocalR=GlobalAffine.TransformPoint(GlobalPosition); 
00052   G4ThreeVector LocalRp=GlobalAffine.TransformAxis(InitMomDir);
00053 
00054   G4ThreeVector itsInitialR=LocalR;
00055   G4ThreeVector itsInitialRp=LocalRp;
00056   
00057   // advance the orbit
00058  
00059   G4ThreeVector itsFinalPoint,itsFinalDir;
00060   
00061   G4ThreeVector yhat(0.,1.,0.);
00062 
00063   G4ThreeVector vhat=LocalRp;
00064   
00065   G4ThreeVector vnorm=vhat.cross(yhat);
00066    
00067   G4double R;
00068   
00069    if(BDSGlobals->GetSynchRescale())
00070     {
00071       G4double B[3];
00072       fPtrMagEqOfMot->GetFieldValue(yIn, B);
00073       R=-(InitMag/GeV)/(0.299792458 * B[1]/tesla) *m;
00074     }
00075   else
00076     {
00077       if(itsBField!=0)
00078         R=-(InitMag/GeV)/(0.299792458 * itsBField/tesla) * m;
00079       else
00080         R=DBL_MAX;
00081     }
00082 
00083  // include the sign of the charge of the particles
00084 
00085   if(  fPtrMagEqOfMot->FCof()<0) R*=-1.;
00086   else if ( fPtrMagEqOfMot->FCof()==0) R=DBL_MAX;
00087 
00088   
00089       G4double Theta   = h/R;
00090 
00091       G4double CosT_ov_2, SinT_ov_2, CosT, SinT;
00092       CosT_ov_2=cos(Theta/2);
00093       SinT_ov_2=sin(Theta/2);
00094 
00095       CosT=(CosT_ov_2*CosT_ov_2)- (SinT_ov_2*SinT_ov_2);
00096       SinT=2*CosT_ov_2*SinT_ov_2;
00097 
00098       BDSLocalRadiusOfCurvature=R;
00099 
00100       itsDist=fabs(R)*(1.-CosT_ov_2);
00101 
00102       G4ThreeVector dPos=R*(SinT*vhat + (1-CosT)*vnorm);
00103         
00104       itsFinalPoint=LocalR+dPos;
00105       itsFinalDir=CosT*vhat +SinT*vnorm;
00106 
00107 
00108   G4ThreeVector GlobalTangent;
00109 
00110   GlobalPosition=LocalAffine.TransformPoint(itsFinalPoint); 
00111   GlobalTangent=LocalAffine.TransformAxis(itsFinalDir);
00112 
00113   GlobalTangent*=InitMag;
00114 
00115   yOut[0]   = GlobalPosition.x(); 
00116   yOut[1]   = GlobalPosition.y(); 
00117   yOut[2]   = GlobalPosition.z(); 
00118   
00119   yOut[3] = GlobalTangent.x();
00120   yOut[4] = GlobalTangent.y();
00121   yOut[5] = GlobalTangent.z();
00122 
00123   
00124 
00125   G4double kappa= - fPtrMagEqOfMot->FCof()* ( itsBGrad) /InitMag; // was ist das? 
00126   if(fabs(kappa)<1.e-12) return; // no gradient
00127 
00128   G4double x1,x1p,y1,y1p,z1p;
00129   //G4double z1;
00130   
00131   G4double NomEnergy = BDSGlobals->GetBeamTotalEnergy();
00132   G4double NomR = -(NomEnergy/GeV)/(0.299792458 * itsBField/tesla) * m;
00133 
00134   G4double NominalPath = sqrt(NomR*NomR - LocalR.z()*LocalR.z()) - fabs(NomR)*cos(itsAngle/2);
00135   
00136   G4double EndNomPath = sqrt(NomR*NomR - itsFinalPoint.z()*itsFinalPoint.z()) - fabs(NomR)*cos(itsAngle/2);
00137 
00138   if(R<0)
00139     {
00140       NominalPath*=-1;
00141       EndNomPath*=-1;
00142     }
00143 
00144   G4double x0=LocalR.x() - NominalPath;
00145   G4double y0=LocalR.y();
00146   G4double z0=LocalR.z();
00147 
00148   G4double theta_in = asin(z0/NomR);
00149   
00150   LocalRp.rotateY(-theta_in);
00151 
00152   G4double xp=LocalRp.x();
00153   G4double yp=LocalRp.y();
00154   G4double zp=LocalRp.z();
00155 
00156   // Save for Synchrotron Radiation calculations:
00157   BDSLocalRadiusOfCurvature=R;
00158   
00159   G4double rootK=sqrt(fabs(kappa*zp));
00160   G4double rootKh=rootK*h*zp;
00161   G4double X11,X12,X21,X22;
00162   G4double Y11,Y12,Y21,Y22;
00163 
00164   if (kappa>0)
00165     {
00166       X11= cos(rootKh);
00167       X12= sin(rootKh)/rootK;
00168       X21=-fabs(kappa)*X12;
00169       X22= X11;
00170       
00171       Y11= cosh(rootKh);
00172       Y12= sinh(rootKh)/rootK;
00173       Y21= fabs(kappa)*Y12;
00174       Y22= Y11;
00175     }
00176   else if (kappa<0)
00177     {
00178       X11= cosh(rootKh);
00179       X12= sinh(rootKh)/rootK;
00180       X21= fabs(kappa)*X12;
00181       X22= X11;
00182       
00183       Y11= cos(rootKh);
00184       Y12= sin(rootKh)/rootK;
00185       Y21= -fabs(kappa)*Y12;
00186       Y22= Y11;
00187     }
00188   else
00189     {
00190       X11 = 1;
00191       X12 = 0;
00192       X21 = 0;
00193       X22 = 1;
00194       Y11 = 1;
00195       Y12 = 0;
00196       Y21 = 0;
00197       Y22 = 1;
00198     }
00199 
00200   x1      = X11*x0 + X12*xp;    
00201   x1p= X21*x0 + X22*xp;
00202 
00203   y1      = Y11*y0 + Y12*yp;    
00204   y1p= Y21*y0 + Y22*yp;
00205 
00206   z1p=sqrt(1 - x1p*x1p -y1p*y1p);
00207   /* 
00208   x1 -=(kappa/ (24*R) ) * h2*h2;
00209   x1p-=(kappa/ (6*R) ) * h*h2;
00210   */
00211   G4double dx=x1-x0;
00212   G4double dy=y1-y0;
00213   // Linear chord length
00214   
00215   LocalR.setX(dx +itsInitialR.x() + EndNomPath - NominalPath);
00216   LocalR.setY(dy + itsInitialR.y());
00217   LocalR.setZ(itsFinalPoint.z());
00218   
00219 
00220   LocalRp.setX(x1p);
00221   LocalRp.setY(y1p);
00222   LocalRp.setZ(z1p);
00223   LocalRp.rotateY(theta_in);
00224   
00225   GlobalPosition=LocalAffine.TransformPoint(LocalR); 
00226   
00227   LocalRp.rotateY(-h/R);
00228   GlobalTangent=LocalAffine.TransformAxis(LocalRp);
00229 
00230 
00231   GlobalTangent*=InitMag;
00232   
00233   // gab: replace += with =
00234   yOut[0] = GlobalPosition.x(); 
00235   yOut[1] = GlobalPosition.y(); 
00236   yOut[2] = GlobalPosition.z(); 
00237   
00238   yOut[3] = GlobalTangent.x();
00239   yOut[4] = GlobalTangent.y();
00240   yOut[5] = GlobalTangent.z();
00241 
00242 }    
00243 
00244 
00245 void myQuadStepper::Stepper( const G4double yInput[],
00246                              const G4double dydx[],
00247                              const G4double hstep,
00248                              G4double yOut[],
00249                              G4double yErr[]      )
00250 {  
00251   const G4int nvar = 6 ;
00252 
00253   for(G4int i=0;i<nvar;i++) yErr[i]=0;
00254 
00255   AdvanceHelix(yInput,0,hstep,yOut);
00256 
00257   return ;
00258 }
00259 
00260 G4double myQuadStepper::DistChord()   const 
00261 {
00262 
00263 return itsDist;
00264   // This is a class method that gives distance of Mid 
00265   //  from the Chord between the Initial and Final points.
00266 }
00267 
00268 myQuadStepper::~myQuadStepper()
00269 {}

Generated on Wed Mar 5 17:25:23 2008 for BDSIM by  doxygen 1.5.3