src/BDSSextStepper.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 // This code implementation is the intellectual property of
00008 // the GEANT4 collaboration.
00009 //
00010 // By copying, distributing or modifying the Program (or any work
00011 // based on the Program) you indicate your acceptance of this statement,
00012 // and all its terms.
00013 //
00014 // $Id: BDSSextStepper.cc,v 1.5 2007/11/14 12:57:06 malton Exp $
00015 // GEANT4 tag $Name:  $
00016 //
00017 #include "BDSSextStepper.hh"
00018 #include "G4ThreeVector.hh"
00019 #include "G4LineSection.hh"
00020 #include "G4TransportationManager.hh"
00021 
00022 extern G4double BDSLocalRadiusOfCurvature;
00023 
00024 BDSSextStepper::BDSSextStepper(G4Mag_EqRhs *EqRhs)
00025   : G4MagIntegratorStepper(EqRhs,6),  // integrate over 6 variables only !!
00026                                       // position & velocity
00027     itsBDblPrime(0.0), itsDist(0.0)
00028 {
00029   fPtrMagEqOfMot = EqRhs;
00030 }
00031 
00032 
00033 void BDSSextStepper::AdvanceHelix( const G4double  yIn[],
00034                                    G4ThreeVector,
00035                                    G4double  h,
00036                                    G4double  ySext[])
00037 {
00038 
00039   const G4double *pIn = yIn+3;
00040   G4ThreeVector v0= G4ThreeVector( pIn[0], pIn[1], pIn[2]);  
00041   G4ThreeVector InitMomDir=v0.unit();
00042 
00043   G4ThreeVector GlobalPosition= G4ThreeVector( yIn[0], yIn[1], yIn[2]);  
00044 
00045   G4double InitMag=v0.mag();
00046 
00047 
00048   G4double kappa=  (-fPtrMagEqOfMot->FCof()*itsBDblPrime) /InitMag;
00049    
00050 
00051 #ifdef DEBUG
00052   G4cout<<"sextupole stepper:"<<G4endl; 
00053   G4cout << "kappa: " << kappa << G4endl;
00054   G4cout << "InitMag: " << InitMag << G4endl;
00055   G4cout << "g'': " <<itsBDblPrime<< G4endl;
00056   G4cout << "fPtrMagEqOfMot->FCof(): " << fPtrMagEqOfMot->FCof() << G4endl << G4endl;
00057   G4cout << "h=: " <<h<< G4endl;
00058 #endif
00059 
00060    if(fabs(kappa)<1.e-12)
00061      {
00062        G4ThreeVector positionMove  = (h/InitMag) * v0;
00063        
00064        ySext[0]   = yIn[0] + positionMove.x(); 
00065        ySext[1]   = yIn[1] + positionMove.y(); 
00066        ySext[2]   = yIn[2] + positionMove.z(); 
00067                                 
00068        ySext[3] = v0.x();
00069        ySext[4] = v0.y();
00070        ySext[5] = v0.z();
00071 
00072        itsDist=0;
00073      }
00074    else 
00075      {      
00076        G4Navigator* SextNavigator=
00077          G4TransportationManager::GetTransportationManager()->
00078          GetNavigatorForTracking();
00079        
00080        G4AffineTransform LocalAffine=SextNavigator->
00081          GetLocalToGlobalTransform();
00082        
00083        G4AffineTransform GlobalAffine=SextNavigator->
00084          GetGlobalToLocalTransform();
00085        G4ThreeVector LocalR=GlobalAffine.TransformPoint(GlobalPosition); 
00086        G4ThreeVector LocalRp=GlobalAffine.TransformAxis(InitMomDir);
00087        
00088        
00089        G4double x0=LocalR.x(); 
00090        G4double y0=LocalR.y();
00091 
00092        // Evaluate field at the approximate midpoint of the step.
00093        x0=x0+LocalRp.x()*h/2;
00094        y0=y0+LocalRp.y()*h/2;
00095        
00096        G4double x02My02=(x0*x0-y0*y0);
00097 
00098        G4double xp=LocalRp.x();
00099        G4double yp=LocalRp.y();
00100        G4double zp=LocalRp.z();
00101 
00102        // local r'' (for curvature)
00103        G4ThreeVector LocalRpp;
00104        LocalRpp.setX(- zp*x02My02);
00105        LocalRpp.setY(2*zp*x0*y0);
00106        LocalRpp.setZ(xp*x02My02-2*yp*x0*y0);
00107 
00108        //G4cout << "LocalRpp: " <<LocalRpp<< G4endl;
00109 
00110        LocalRpp*=kappa/2; // 2 is actually a 2! factor.
00111        // determine effective curvature
00112        G4double R_1 = LocalRpp.mag();
00113 
00114 
00115        if(R_1>0.)
00116          {    
00117            G4double h2=h*h;
00118            // chord distance (simple quadratic approx)
00119            itsDist= h2*R_1/8;
00120 
00121            // Save for Synchrotron Radiation calculations:
00122            BDSLocalRadiusOfCurvature=1./R_1;
00123            
00124            G4double dx=LocalRp.x()*h + LocalRpp.x()*h2 /2.; 
00125            G4double dy=LocalRp.y()*h + LocalRpp.y()*h2 /2.;
00126 
00127            G4double dz=sqrt(h2*(1.-h2*R_1*R_1/12)-dx*dx-dy*dy);
00128            // check for precision problems
00129            G4double ScaleFac=(dx*dx+dy*dy+dz*dz)/h2;
00130            if(ScaleFac>1.0000001)
00131              {
00132                ScaleFac=sqrt(ScaleFac);
00133                dx/=ScaleFac;
00134                dy/=ScaleFac;
00135                dz/=ScaleFac;
00136              }
00137 
00138          
00139            LocalR.setX(LocalR.x()+dx);
00140            LocalR.setY(LocalR.y()+dy);
00141            LocalR.setZ(LocalR.z()+dz);
00142           
00143            LocalRp = LocalRp+ h*LocalRpp;
00144 
00145          }
00146        else
00147          {LocalR += h*LocalRp;}
00148        
00149        GlobalPosition=LocalAffine.TransformPoint(LocalR); 
00150        G4ThreeVector GlobalTangent=LocalAffine.TransformAxis(LocalRp)*InitMag;
00151        
00152        ySext[0]   = GlobalPosition.x(); 
00153        ySext[1]   = GlobalPosition.y(); 
00154        ySext[2]   = GlobalPosition.z(); 
00155                                 
00156        ySext[3] = GlobalTangent.x();
00157        ySext[4] = GlobalTangent.y();
00158        ySext[5] = GlobalTangent.z();
00159      }
00160 }
00161 
00162 void BDSSextStepper::Stepper( const G4double yInput[],
00163                              const G4double[],
00164                              const G4double hstep,
00165                              G4double yOut[],
00166                              G4double yErr[]      )
00167 {  
00168   const G4int nvar = 6 ;
00169   
00170   G4int i;
00171  
00172   G4double yTemp[7], yIn[7];
00173   
00174   //  Saving yInput because yInput and yOut can be aliases for same array
00175   
00176   for(i=0;i<nvar;i++) yIn[i]=yInput[i];
00177   
00178   G4double h = hstep * 0.5; 
00179   
00180   // Do two half steps
00181   AdvanceHelix(yIn,   (G4ThreeVector)0,  h, yTemp);
00182   AdvanceHelix(yTemp, (G4ThreeVector)0, h, yOut); 
00183   
00184   // Do a full Step
00185   h = hstep ;
00186   AdvanceHelix(yIn, (G4ThreeVector)0, h, yTemp); 
00187   
00188   for(i=0;i<nvar;i++) yErr[i] = yOut[i] - yTemp[i] ;
00189 }
00190 
00191 G4double BDSSextStepper::DistChord()   const 
00192 {
00193   return itsDist;
00194 }
00195 
00196 BDSSextStepper::~BDSSextStepper()
00197 {}

Generated on 27 Aug 2013 for BDSIM by  doxygen 1.4.7