G4RKG3_Stepper.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id: G4RKG3_Stepper.cc 69786 2013-05-15 09:38:51Z gcosmo $
00028 //
00029 // -------------------------------------------------------------------
00030 
00031 #include "G4RKG3_Stepper.hh"
00032 #include "G4LineSection.hh"
00033 #include "G4Mag_EqRhs.hh"
00034 
00035 G4RKG3_Stepper::G4RKG3_Stepper(G4Mag_EqRhs *EqRhs)
00036   : G4MagIntegratorStepper(EqRhs,6), hStep(0.)
00037 {
00038 }
00039 
00040 G4RKG3_Stepper::~G4RKG3_Stepper()
00041 {
00042 }
00043 
00044 void G4RKG3_Stepper::Stepper(  const G4double yInput[8],
00045                                const G4double dydx[6],
00046                                      G4double Step,
00047                                      G4double yOut[8],
00048                                      G4double yErr[])
00049 {
00050    G4double  B[3];
00051    G4int nvar = 6 ;
00052    G4int i;
00053    G4double  by15 = 1. / 15. ; // was  0.066666666 ;
00054 
00055    G4double yTemp[8], dydxTemp[6], yIn[8] ;
00056    //  Saving yInput because yInput and yOut can be aliases for same array
00057    for(i=0;i<nvar;i++) yIn[i]=yInput[i];
00058    yIn[6] = yInput[6];
00059    yIn[7] = yInput[7];
00060    G4double h = Step * 0.5; 
00061    hStep=Step;
00062    // Do two half steps
00063 
00064    StepNoErr(yIn, dydx,h, yTemp,B) ;
00065    
00066    //Store Bfld for DistChord Calculation
00067    for(i=0;i<3;i++)BfldIn[i]=B[i];
00068 
00069    //   RightHandSide(yTemp,dydxTemp) ;
00070 
00071    GetEquationOfMotion()->EvaluateRhsGivenB(yTemp,B,dydxTemp) ;  
00072    StepNoErr(yTemp,dydxTemp,h,yOut,B);      
00073         
00074    // Store midpoint, chord calculation
00075                                  
00076    fyMidPoint = G4ThreeVector( yTemp[0],  yTemp[1],  yTemp[2]); 
00077 
00078    // Do a full Step
00079 
00080    h *= 2 ;
00081    StepNoErr(yIn,dydx,h,yTemp,B); 
00082    for(i=0;i<nvar;i++)
00083    {
00084       yErr[i] = yOut[i] - yTemp[i] ;
00085       yOut[i] += yErr[i]*by15 ;          // Provides 5th order of accuracy
00086    }
00087 
00088    //Store values for DistChord method
00089 
00090    fyInitial = G4ThreeVector( yIn[0],   yIn[1],   yIn[2]);
00091    fpInitial = G4ThreeVector( yIn[3],   yIn[4],   yIn[5]);
00092    fyFinal   = G4ThreeVector( yOut[0],  yOut[1],  yOut[2]); 
00093   
00094    // NormaliseTangentVector( yOut );  // Deleted
00095 }
00096 
00097 // ---------------------------------------------------------------------------
00098 
00099 // Integrator for RK from G3 with evaluation of error in solution and delta
00100 // geometry based on naive similarity with the case of uniform magnetic field.
00101 // B1[3] is input  and is the first magnetic field values
00102 // B2[3] is output and is the final magnetic field values.
00103 
00104 void G4RKG3_Stepper::StepWithEst( const G4double*,
00105                                   const G4double*,
00106                                         G4double,
00107                                         G4double*,
00108                                         G4double&,
00109                                         G4double&,
00110                                   const G4double*,
00111                                         G4double* )
00112    
00113 {
00114   G4Exception("G4RKG3_Stepper::StepWithEst()", "GeomField0001",
00115               FatalException, "Method no longer used.");
00116 }
00117 
00118 // -----------------------------------------------------------------
00119 
00120 
00121 // Integrator RK Stepper from G3 with only two field evaluation per Step. 
00122 // It is used in propagation initial Step by small substeps after solution 
00123 // error and delta geometry considerations. B[3] is magnetic field which 
00124 // is passed from substep to substep.
00125 
00126 void G4RKG3_Stepper::StepNoErr(const G4double tIn[8],
00127                                const G4double dydx[6],
00128                                      G4double Step,
00129                                      G4double tOut[8],
00130                                      G4double B[3]      )     // const
00131    
00132 { 
00133   
00134   //  Copy and edit the routine above, to delete alpha2, beta2, ...
00135    G4double K1[7],K2[7],K3[7],K4[7] ;
00136    G4double tTemp[8], yderiv[6] ;
00137 
00138   // Need Momentum value to give correct values to the coefficients in equation
00139   // Integration on unit velocity, but  tIn[3,4,5] is momentum 
00140    G4double mom,inverse_mom;
00141    G4int i ;
00142    const G4double c1=0.5,c2=0.125,c3=1./6.;
00143   
00144    // GetEquationOfMotion()->EvaluateRhsReturnB(tIn,dydx,B1) ;
00145    // Correction for momentum not a velocity
00146    // Need the protection !!! must be not zero 
00147      mom=std::sqrt(tIn[3]*tIn[3]+tIn[4]*tIn[4]+tIn[5]*tIn[5]); 
00148      inverse_mom=1./mom;    
00149    for(i=0;i<3;i++)
00150    {
00151       K1[i] = Step * dydx[i+3]*inverse_mom;
00152       tTemp[i] = tIn[i] + Step*(c1*tIn[i+3]*inverse_mom + c2*K1[i]) ;
00153       tTemp[i+3] = tIn[i+3] + c1*K1[i]*mom ;
00154      
00155    }
00156     
00157    GetEquationOfMotion()->EvaluateRhsReturnB(tTemp,yderiv,B) ;
00158     
00159       
00160    for(i=0;i<3;i++)
00161    {
00162       K2[i] = Step * yderiv[i+3]*inverse_mom;
00163       tTemp[i+3] = tIn[i+3] + c1*K2[i]*mom ;
00164    }
00165    
00166    //  Given B, calculate yderiv !
00167    GetEquationOfMotion()->EvaluateRhsGivenB(tTemp,B,yderiv) ;  
00168  
00169    for(i=0;i<3;i++)
00170    {
00171       K3[i] = Step * yderiv[i+3]*inverse_mom;
00172       tTemp[i] = tIn[i] + Step*(tIn[i+3]*inverse_mom + c1*K3[i]) ;
00173       tTemp[i+3] = tIn[i+3] + K3[i]*mom ;
00174    }
00175    
00176 
00177    //  Calculates y-deriv(atives) & returns B too!
00178    GetEquationOfMotion()->EvaluateRhsReturnB(tTemp,yderiv,B) ;  
00179     
00180 
00181    for(i=0;i<3;i++)        // Output trajectory vector
00182    {
00183       K4[i] = Step * yderiv[i+3]*inverse_mom;
00184       tOut[i] = tIn[i] + Step*(tIn[i+3]*inverse_mom+ (K1[i] + K2[i] + K3[i])*c3) ;
00185       tOut[i+3] = tIn[i+3] + mom*(K1[i] + 2*K2[i] + 2*K3[i] +K4[i])*c3 ;
00186    }
00187    tOut[6] = tIn[6];
00188    tOut[7] = tIn[7];
00189    // NormaliseTangentVector( tOut );
00190   
00191 
00192 }
00193 
00194 
00195 // ---------------------------------------------------------------------------
00196  
00197  G4double G4RKG3_Stepper::DistChord()   const 
00198  {
00199    // Soon: must check whether h/R > 2 pi  !!
00200    //  Method below is good only for < 2 pi
00201    G4double distChord,distLine;
00202    
00203    if (fyInitial != fyFinal) {
00204       distLine= G4LineSection::Distline(fyMidPoint,fyInitial,fyFinal );
00205   
00206         distChord = distLine;
00207    }else{
00208       distChord = (fyMidPoint-fyInitial).mag();
00209    }
00210  
00211   
00212    return distChord;
00213    
00214  }
00215 

Generated on Mon May 27 17:49:44 2013 for Geant4 by  doxygen 1.4.7