/scratch0/jsnuveri/BDSIM/BDSIMgit/bdsim/src/BDSOctStepper.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 // This code implementation is the intellectual property of
00007 // the GEANT4 collaboration.
00008 //
00009 // By copying, distributing or modifying the Program (or any work
00010 // based on the Program) you indicate your acceptance of this statement,
00011 // and all its terms.
00012 //
00013 // $Id: BDSOctStepper.cc,v 1.6 2007/05/10 16:23:22 malton Exp $
00014 // GEANT4 tag $Name:  $
00015 //
00016 #include "BDSOctStepper.hh"
00017 #include "G4ThreeVector.hh"
00018 #include "G4TransportationManager.hh"
00019 
00020 extern G4double BDSLocalRadiusOfCurvature;
00021 
00022 BDSOctStepper::BDSOctStepper(G4Mag_EqRhs *EqRhs)
00023   : G4MagIntegratorStepper(EqRhs,6),  // integrate over 6 variables only !!
00024                                       // position & velocity
00025     itsBTrpPrime(0.0), itsDist(0.0)
00026 {
00027   fPtrMagEqOfMot = EqRhs;
00028 }
00029 
00030 
00031 void BDSOctStepper::AdvanceHelix( const G4double  yIn[],
00032                                    G4ThreeVector,
00033                                    G4double  h,
00034                                    G4double  yOct[])
00035 {
00036   const G4double *pIn = yIn+3;
00037   G4ThreeVector v0= G4ThreeVector( pIn[0], pIn[1], pIn[2]);  
00038   G4ThreeVector InitMomDir=v0.unit();
00039 
00040   G4ThreeVector GlobalPosition= G4ThreeVector( yIn[0], yIn[1], yIn[2]);  
00041   G4double InitMag=v0.mag();
00042   G4double kappa=  -fPtrMagEqOfMot->FCof()*itsBTrpPrime/InitMag;
00043 
00044   // relevant momentum scale is p_z, not P_tot:
00045   // check that the approximations are valid, else do a linear step:
00046   if(fabs(kappa)<1.e-20)
00047     {
00048       G4ThreeVector positionMove  = (h/InitMag) * v0;
00049       
00050       yOct[0]   = yIn[0] + positionMove.x(); 
00051       yOct[1]   = yIn[1] + positionMove.y(); 
00052       yOct[2]   = yIn[2] + positionMove.z(); 
00053                                 
00054       yOct[3] = v0.x();
00055       yOct[4] = v0.y();
00056       yOct[5] = v0.z();
00057 
00058       itsDist=0;
00059     }
00060   else 
00061     {      
00062       G4Navigator* OctNavigator=
00063         G4TransportationManager::GetTransportationManager()->
00064         GetNavigatorForTracking();
00065 
00066       G4AffineTransform LocalAffine=OctNavigator->GetLocalToGlobalTransform();
00067 
00068 
00069       // gab_dec03>>
00070       // position 
00071       //G4ThreeVector LocalR = OctNavigator->GetCurrentLocalCoordinate();
00072       // position derivative r' (normalised to unity)
00073       //G4ThreeVector LocalRp= (OctNavigator->ComputeLocalAxis(v0)).unit();
00074 
00075       G4AffineTransform GlobalAffine=OctNavigator->GetGlobalToLocalTransform();
00076       G4ThreeVector LocalR=GlobalAffine.TransformPoint(GlobalPosition); 
00077       G4ThreeVector LocalRp=GlobalAffine.TransformAxis(InitMomDir);
00078       // gab_dec03<<
00079 
00080 
00081 
00082       G4double x0=LocalR.x(); 
00083       G4double y0=LocalR.y();
00084       G4double z0=LocalR.z();
00085 
00086       //G4double x02My02=(x0*x0-y0*y0);
00087 
00088       G4double xp=LocalRp.x();
00089       G4double yp=LocalRp.y();
00090       G4double zp=LocalRp.z();
00091 
00092       G4double y3fac=y0*(y0*y0-3*x0*x0);
00093       G4double x3fac=x0*(x0*x0-3*y0*y0);
00094       
00095       // local r'' (for curvature)
00096       G4ThreeVector LocalRpp;
00097       // extra minus signs were because x,y_machine = - x_,-y_geant_world
00098       // New CVS version of BDSIM uses +x, +y in geant world
00099       /*
00100       LocalRpp.setX(zp*x3fac);
00101       LocalRpp.setY(zp*y3fac);
00102       LocalRpp.setZ(- xp*x3fac - yp*y3fac);
00103       */
00104       LocalRpp.setX(-zp*x3fac);
00105       LocalRpp.setY(-zp*y3fac);
00106       LocalRpp.setZ( xp*x3fac + yp*y3fac);
00107 
00108       LocalRpp*=kappa/6; // 6 is actually a 3! factor.;
00109 
00110       // determine effective curvature
00111       G4double R_1 = LocalRpp.mag();
00112       if(R_1>0.)
00113         {
00114           // Save for Synchrotron Radiation calculations:
00115           BDSLocalRadiusOfCurvature=1/R_1;
00116 
00117           // chord distance (simple quadratic approx)
00118           G4double h2=h*h;
00119           itsDist= h2*R_1/8;
00120 
00121           G4double dx=LocalRp.x()*h + LocalRpp.x()*h2/2; 
00122           G4double dy=LocalRp.y()*h + LocalRpp.y()*h2/2; 
00123 
00124           G4double dz=sqrt(h2*(1.-h2*R_1*R_1/12)-dx*dx-dy*dy);
00125           // check for precision problems
00126           G4double ScaleFac=(dx*dx+dy*dy+dz*dz)/h2;
00127           if(ScaleFac>1.0000001)
00128             {
00129               ScaleFac=sqrt(ScaleFac);
00130               dx/=ScaleFac;
00131               dy/=ScaleFac;
00132               dz/=ScaleFac;
00133             }
00134 
00135           LocalR.setX(x0+dx);
00136           LocalR.setY(y0+dy);
00137           LocalR.setZ(z0+dz);
00138 
00139           LocalRp = LocalRp+ h*LocalRpp;
00140         }
00141       else
00142         {LocalR += h*LocalRp;}
00143        
00144       GlobalPosition=LocalAffine.TransformPoint(LocalR); 
00145       G4ThreeVector GlobalTangent=LocalAffine.TransformAxis(LocalRp)*InitMag;
00146       
00147       yOct[0]   = GlobalPosition.x(); 
00148       yOct[1]   = GlobalPosition.y(); 
00149       yOct[2]   = GlobalPosition.z(); 
00150                                 
00151       yOct[3] = GlobalTangent.x();
00152       yOct[4] = GlobalTangent.y();
00153       yOct[5] = GlobalTangent.z();
00154     }
00155 }
00156 
00157 void BDSOctStepper::Stepper( const G4double yInput[],
00158                             const G4double[],
00159                             const G4double hstep,
00160                             G4double yOut[],
00161                             G4double yErr[]      )
00162 {  
00163   G4int i;
00164   const G4int nvar = 6 ;
00165   
00166   //const G4double *pIn = yInput+3;
00167   //G4ThreeVector v0= G4ThreeVector( pIn[0], pIn[1], pIn[2]);  
00168   //G4double InitMag=v0.mag();
00169   //G4double kappa=  -fPtrMagEqOfMot->FCof()*itsBTrpPrime/InitMag;
00170   
00171   G4double yTemp[7], yIn[7];
00172   
00173   //  Saving yInput because yInput and yOut can be aliases for same array
00174   
00175   for(i=0;i<nvar;i++) yIn[i]=yInput[i];
00176   
00177   G4double h = hstep * 0.5; 
00178   
00179   // Do two half steps
00180   AdvanceHelix(yIn,   (G4ThreeVector)0,  h, yTemp);
00181   AdvanceHelix(yTemp, (G4ThreeVector)0, h, yOut); 
00182   
00183   // Do a full Step
00184   h = hstep ;
00185   AdvanceHelix(yIn, (G4ThreeVector)0, h, yTemp); 
00186   
00187   for(i=0;i<nvar;i++) {
00188     yErr[i] = yOut[i] - yTemp[i] ;
00189     // if error small, set error to 0
00190     // this is done to prevent Geant4 going to smaller and smaller steps
00191     // ideally use some of the global constants instead of hardcoding here
00192     // could look at step size as well instead.
00193     if (std::abs(yErr[i]) < 1e-7) yErr[i] = 0;
00194   }
00195 }
00196 
00197 G4double BDSOctStepper::DistChord()   const 
00198 {
00199   return itsDist;
00200 }
00201 
00202 BDSOctStepper::~BDSOctStepper()
00203 {}

Generated on 28 Jun 2015 for BDSIM by  doxygen 1.4.7