src/BDSDecStepper.cc

00001 /* BDSIM code.    Version 1.0
00002    Author: Grahame A. Blair, Royal Holloway, Univ. of London.
00003    Last modified 25.12.2003
00004    Copyright (c) 2003 by G.A.Blair.  ALL RIGHTS RESERVED. 
00005 */
00006 
00007 #include "BDSDecStepper.hh"
00008 #include "G4ThreeVector.hh"
00009 #include "G4LineSection.hh"
00010 #include "G4TransportationManager.hh"
00011 
00012 extern G4double BDSLocalRadiusOfCurvature;
00013 
00014 BDSDecStepper::BDSDecStepper(G4Mag_EqRhs *EqRhs)
00015   : G4MagIntegratorStepper(EqRhs,6),   // integrate over 6 variables only !!
00016                                        // position & velocity
00017     itsBQuadPrime(0.0), itsDist(0.0)
00018 {
00019   fPtrMagEqOfMot = EqRhs;
00020 }
00021 
00022 
00023 void BDSDecStepper::AdvanceHelix( const G4double  yIn[],
00024                                    G4ThreeVector,
00025                                    G4double  h,
00026                                    G4double  yDec[])
00027 {
00028   const G4double *pIn = yIn+3;
00029   G4ThreeVector v0= G4ThreeVector( pIn[0], pIn[1], pIn[2]);  
00030   G4ThreeVector InitMomDir=v0.unit();
00031 
00032   G4ThreeVector GlobalPosition= G4ThreeVector( yIn[0], yIn[1], yIn[2]);  
00033   G4double InitMag=v0.mag();
00034   G4double kappa=  -fPtrMagEqOfMot->FCof()*itsBQuadPrime/InitMag;
00035 
00036   // relevant momentum scale is p_z, not P_tot:
00037   // check that the approximations are valid, else do a linear step:
00038   if(fabs(kappa)<1.e-20)
00039     {
00040       G4ThreeVector positionMove  = (h/InitMag) * v0;
00041       
00042       yDec[0]   = yIn[0] + positionMove.x(); 
00043       yDec[1]   = yIn[1] + positionMove.y(); 
00044       yDec[2]   = yIn[2] + positionMove.z(); 
00045                                 
00046       yDec[3] = v0.x();
00047       yDec[4] = v0.y();
00048       yDec[5] = v0.z();
00049 
00050       itsDist=0;
00051     }
00052   else 
00053     {      
00054       G4Navigator* DecNavigator=
00055         G4TransportationManager::GetTransportationManager()->
00056         GetNavigatorForTracking();
00057 
00058       G4AffineTransform LocalAffine=DecNavigator->GetLocalToGlobalTransform();
00059 
00060 
00061       // gab_dec03>>
00062       // position 
00063       // gab_dec03
00064       //      G4ThreeVector LocalR = DecNavigator->GetCurrentLocalCoordinate();
00065       // position derivative r' (normalised to unity)
00066       //      G4ThreeVector LocalRp= (DecNavigator->ComputeLocalAxis(v0)).unit();
00067 
00068       G4AffineTransform GlobalAffine=DecNavigator->GetGlobalToLocalTransform();
00069       G4ThreeVector LocalR=GlobalAffine.TransformPoint(GlobalPosition); 
00070       G4ThreeVector LocalRp=GlobalAffine.TransformAxis(InitMomDir);
00071       // gab_dec03<<
00072 
00073 
00074       G4double x0=LocalR.x(); 
00075       G4double y0=LocalR.y();
00076     
00077       G4double x02My02=(x0*x0-y0*y0);
00078 
00079       G4double xp=LocalRp.x();
00080       G4double yp=LocalRp.y();
00081       G4double zp=LocalRp.z();
00082 
00083       G4double Bx=4.0*x0*y0*(x02My02);
00084       G4double By=pow(x0,4)-6.0*x0*x0*y0*y0+pow(y0,4);
00085 
00086       // local r'' (for curvature)
00087       G4ThreeVector LocalRpp;
00088       // extra minus signs because x,y_machine = - x_,-y_geant_world
00089       LocalRpp.setX(zp*By);
00090       LocalRpp.setY(-zp*Bx);
00091       LocalRpp.setZ( xp*By - yp*Bx);
00092       
00093       LocalRpp*=kappa/24; // 24 is actually a 4! factor.;
00094 
00095       // determine effective curvature
00096       G4double R_1 = LocalRpp.mag();
00097       if(R_1>0.)
00098         {
00099           // Save for Synchrotron Radiation calculations:
00100           BDSLocalRadiusOfCurvature=1/R_1;
00101 
00102           // chord distance (simple quadratic approx)
00103           G4double h2=h*h;
00104           itsDist= h2*R_1/8;
00105 
00106           G4double dx=LocalRp.x()*h + LocalRpp.x()*h2/2; 
00107           G4double dy=LocalRp.y()*h + LocalRpp.y()*h2/2; 
00108 
00109           G4double dz=sqrt(h2*(1.-h2*R_1*R_1/12)-dx*dx-dy*dy);
00110           // check for precision problems
00111           G4double ScaleFac=(dx*dx+dy*dy+dz*dz)/h2;
00112           if(ScaleFac>1.0000001)
00113             {
00114               ScaleFac=sqrt(ScaleFac);
00115               dx/=ScaleFac;
00116               dy/=ScaleFac;
00117               dz/=ScaleFac;
00118             }
00119 
00120           LocalR.setX(LocalR.x()+dx);
00121           LocalR.setY(LocalR.y()+dy);
00122           LocalR.setZ(LocalR.z()+dz);
00123 
00124           LocalRp = LocalRp+ h*LocalRpp;
00125         }
00126       else
00127         {LocalR += h*LocalRp;}
00128        
00129       GlobalPosition=LocalAffine.TransformPoint(LocalR); 
00130       G4ThreeVector GlobalTangent=LocalAffine.TransformAxis(LocalRp)*InitMag;
00131       
00132       yDec[0]   = GlobalPosition.x(); 
00133       yDec[1]   = GlobalPosition.y(); 
00134       yDec[2]   = GlobalPosition.z(); 
00135                                 
00136       yDec[3] = GlobalTangent.x();
00137       yDec[4] = GlobalTangent.y();
00138       yDec[5] = GlobalTangent.z();
00139     }
00140 }
00141 
00142 void BDSDecStepper::Stepper( const G4double yInput[],
00143                             const G4double[],
00144                             const G4double hstep,
00145                             G4double yOut[],
00146                             G4double yErr[]      )
00147 {  
00148   const G4int nvar = 6 ;
00149   
00150   G4int i;
00151   for(i=0;i<nvar;i++) yErr[i]=0;
00152   AdvanceHelix(yInput,(G4ThreeVector)0,hstep,yOut);
00153 }
00154 
00155 G4double BDSDecStepper::DistChord()   const 
00156 {
00157   return itsDist;
00158 }
00159 
00160 BDSDecStepper::~BDSDecStepper()
00161 {}

Generated on 27 Aug 2013 for BDSIM by  doxygen 1.4.7