G4ConstRK4.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: G4ConstRK4.cc 69786 2013-05-15 09:38:51Z gcosmo $
00028 //
00029 //
00030 // - 18.09.2008 - J.Apostolakis, T.Nikitina - Created
00031 // -------------------------------------------------------------------
00032 
00033 #include "G4ConstRK4.hh"
00034 #include "G4ThreeVector.hh"
00035 #include "G4LineSection.hh"
00036 
00038 //
00039 // Constructor sets the number of *State* variables (default = 8)
00040 //   The number of variables integrated is always 6
00041 
00042 G4ConstRK4::G4ConstRK4(G4Mag_EqRhs* EqRhs, G4int numStateVariables)
00043   : G4MagErrorStepper(EqRhs, 6, numStateVariables)
00044 {
00045   // const G4int numberOfVariables= 6;
00046   if( numStateVariables < 8 ) 
00047   {
00048     std::ostringstream message;
00049     message << "The number of State variables at least 8 " << G4endl
00050             << "Instead it is - numStateVariables= " << numStateVariables;
00051     G4Exception("G4ConstRK4::G4ConstRK4()", "GeomField0002",
00052                 FatalException, message, "Use another Stepper!");
00053   }
00054 
00055   fEq = EqRhs;
00056   yMiddle  = new G4double[8];
00057   dydxMid  = new G4double[8];
00058   yInitial = new G4double[8];
00059   yOneStep = new G4double[8];
00060 
00061   dydxm = new G4double[8];
00062   dydxt = new G4double[8]; 
00063   yt    = new G4double[8]; 
00064   Field[0]=0.; Field[1]=0.; Field[2]=0.;
00065 }
00066 
00068 //
00069 // Destructor
00070 
00071 G4ConstRK4::~G4ConstRK4()
00072 {
00073    delete [] yMiddle;
00074    delete [] dydxMid;
00075    delete [] yInitial;
00076    delete [] yOneStep;
00077    delete [] dydxm;
00078    delete [] dydxt;
00079    delete [] yt;
00080 }
00081 
00083 //
00084 // Given values for the variables y[0,..,n-1] and their derivatives
00085 // dydx[0,...,n-1] known at x, use the classical 4th Runge-Kutta
00086 // method to advance the solution over an interval h and return the
00087 // incremented variables as yout[0,...,n-1], which is not a distinct
00088 // array from y. The user supplies the routine RightHandSide(x,y,dydx),
00089 // which returns derivatives dydx at x. The source is routine rk4 from
00090 // NRC p. 712-713 .
00091 
00092 void G4ConstRK4::DumbStepper( const G4double  yIn[],
00093                               const G4double  dydx[],
00094                                     G4double  h,
00095                                     G4double  yOut[])
00096 {
00097    G4double  hh = h*0.5 , h6 = h/6.0  ;
00098    
00099    // 1st Step K1=h*dydx
00100    yt[5] = yIn[5] + hh*dydx[5] ;
00101    yt[4] = yIn[4] + hh*dydx[4] ;
00102    yt[3] = yIn[3] + hh*dydx[3] ;
00103    yt[2] = yIn[2] + hh*dydx[2] ;
00104    yt[1] = yIn[1] + hh*dydx[1] ;
00105    yt[0] = yIn[0] + hh*dydx[0] ;
00106    RightHandSideConst(yt,dydxt) ;        
00107 
00108    // 2nd Step K2=h*dydxt
00109    yt[5] = yIn[5] + hh*dydxt[5] ;
00110    yt[4] = yIn[4] + hh*dydxt[4] ;
00111    yt[3] = yIn[3] + hh*dydxt[3] ;
00112    yt[2] = yIn[2] + hh*dydxt[2] ;
00113    yt[1] = yIn[1] + hh*dydxt[1] ;
00114    yt[0] = yIn[0] + hh*dydxt[0] ;
00115    RightHandSideConst(yt,dydxm) ;     
00116 
00117    // 3rd Step K3=h*dydxm
00118    // now dydxm=(K2+K3)/h
00119    yt[5] = yIn[5] + h*dydxm[5] ;
00120    dydxm[5] += dydxt[5] ;  
00121    yt[4] = yIn[4] + h*dydxm[4] ;
00122    dydxm[4] += dydxt[4] ;  
00123    yt[3] = yIn[3] + h*dydxm[3] ;
00124    dydxm[3] += dydxt[3] ;  
00125    yt[2] = yIn[2] + h*dydxm[2] ;
00126    dydxm[2] += dydxt[2] ;  
00127    yt[1] = yIn[1] + h*dydxm[1] ;
00128    dydxm[1] += dydxt[1] ;  
00129    yt[0] = yIn[0] + h*dydxm[0] ;
00130    dydxm[0] += dydxt[0] ;  
00131    RightHandSideConst(yt,dydxt) ;   
00132 
00133    // 4th Step K4=h*dydxt
00134    yOut[5] = yIn[5]+h6*(dydx[5]+dydxt[5]+2.0*dydxm[5]);
00135    yOut[4] = yIn[4]+h6*(dydx[4]+dydxt[4]+2.0*dydxm[4]);
00136    yOut[3] = yIn[3]+h6*(dydx[3]+dydxt[3]+2.0*dydxm[3]);
00137    yOut[2] = yIn[2]+h6*(dydx[2]+dydxt[2]+2.0*dydxm[2]);
00138    yOut[1] = yIn[1]+h6*(dydx[1]+dydxt[1]+2.0*dydxm[1]);
00139    yOut[0] = yIn[0]+h6*(dydx[0]+dydxt[0]+2.0*dydxm[0]);
00140    
00141 }  // end of DumbStepper ....................................................
00142 
00144 //
00145 // Stepper
00146 
00147 void
00148 G4ConstRK4::Stepper( const G4double yInput[],
00149                      const G4double dydx[],
00150                            G4double hstep,
00151                            G4double yOutput[],
00152                            G4double yError [] )
00153 {
00154    const G4int nvar = 6;  // number of variables integrated
00155    const G4int maxvar= GetNumberOfStateVariables();
00156 
00157    // Correction for Richardson extrapolation
00158    G4double  correction = 1. / ( (1 << IntegratorOrder()) -1 );
00159 
00160    G4int i;
00161    
00162    // Saving yInput because yInput and yOutput can be aliases for same array
00163    for (i=0;    i<maxvar; i++) { yInitial[i]= yInput[i]; }
00164  
00165    // Must copy the part of the state *not* integrated to the output
00166    for (i=nvar; i<maxvar; i++) { yOutput[i]=  yInput[i]; }
00167 
00168    // yInitial[7]= yInput[7];  //  The time is typically needed
00169    yMiddle[7]  = yInput[7];   // Copy the time from initial value 
00170    yOneStep[7] = yInput[7];   // As it contributes to final value of yOutput ?
00171    // yOutput[7] = yInput[7];  // -> dumb stepper does it too for RK4
00172    yError[7] = 0.0;         
00173 
00174    G4double halfStep = hstep * 0.5; 
00175 
00176    // Do two half steps
00177    //
00178    GetConstField(yInitial,Field);
00179    DumbStepper  (yInitial,  dydx,   halfStep, yMiddle);
00180    RightHandSideConst(yMiddle, dydxMid);    
00181    DumbStepper  (yMiddle, dydxMid, halfStep, yOutput); 
00182 
00183    // Store midpoint, chord calculation
00184    //
00185    fMidPoint = G4ThreeVector( yMiddle[0],  yMiddle[1],  yMiddle[2]); 
00186 
00187    // Do a full Step
00188    //
00189    DumbStepper(yInitial, dydx, hstep, yOneStep);
00190    for(i=0;i<nvar;i++)
00191    {
00192       yError [i] = yOutput[i] - yOneStep[i] ;
00193       yOutput[i] += yError[i]*correction ;
00194         // Provides accuracy increased by 1 order via the 
00195         // Richardson extrapolation  
00196    }
00197 
00198    fInitialPoint = G4ThreeVector( yInitial[0], yInitial[1], yInitial[2]); 
00199    fFinalPoint   = G4ThreeVector( yOutput[0],  yOutput[1],  yOutput[2]); 
00200 
00201    return;
00202 }
00203 
00205 //
00206 // Estimate the maximum distance from the curve to the chord
00207 //
00208 // We estimate this using the distance of the midpoint to chord.
00209 // The method below is good only for angle deviations < 2 pi;
00210 // this restriction should not be a problem for the Runge Kutta methods, 
00211 // which generally cannot integrate accurately for large angle deviations
00212 
00213 G4double G4ConstRK4::DistChord() const 
00214 {
00215   G4double distLine, distChord; 
00216 
00217   if (fInitialPoint != fFinalPoint)
00218   {
00219      distLine= G4LineSection::Distline( fMidPoint, fInitialPoint, fFinalPoint );
00220        // This is a class method that gives distance of Mid 
00221        // from the Chord between the Initial and Final points
00222      distChord = distLine;
00223   }
00224   else
00225   {
00226      distChord = (fMidPoint-fInitialPoint).mag();
00227   }
00228   return distChord;
00229 }

Generated on Mon May 27 17:47:57 2013 for Geant4 by  doxygen 1.4.7