source: trunk/source/geometry/magneticfield/src/G4ConstRK4.cc

Last change on this file was 1340, checked in by garnier, 14 years ago

update ti head

File size: 8.1 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer                                           *
4// *                                                                  *
5// * The  Geant4 software  is  copyright of the Copyright Holders  of *
6// * the Geant4 Collaboration.  It is provided  under  the terms  and *
7// * conditions of the Geant4 Software License,  included in the file *
8// * LICENSE and available at  http://cern.ch/geant4/license .  These *
9// * include a list of copyright holders.                             *
10// *                                                                  *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work  make  any representation or  warranty, express or implied, *
14// * regarding  this  software system or assume any liability for its *
15// * use.  Please see the license in the file  LICENSE  and URL above *
16// * for the full disclaimer and the limitation of liability.         *
17// *                                                                  *
18// * This  code  implementation is the result of  the  scientific and *
19// * technical work of the GEANT4 collaboration.                      *
20// * By using,  copying,  modifying or  distributing the software (or *
21// * any work based  on the software)  you  agree  to acknowledge its *
22// * use  in  resulting  scientific  publications,  and indicate your *
23// * acceptance of all terms of the Geant4 Software license.          *
24// ********************************************************************
25//
26//
27// $Id: G4ConstRK4.cc,v 1.5 2010/09/10 15:51:10 japost Exp $
28// GEANT4 tag $Name: field-V09-03-03 $
29//
30//
31// - 18.09.2008 - J.Apostolakis, T.Nikitina - Created
32// -------------------------------------------------------------------
33
34#include "G4ConstRK4.hh"
35#include "G4ThreeVector.hh"
36#include "G4LineSection.hh"
37
38//////////////////////////////////////////////////////////////////
39//
40// Constructor sets the number of *State* variables (default = 8)
41//   The number of variables integrated is always 6
42
43G4ConstRK4::G4ConstRK4(G4Mag_EqRhs* EqRhs, G4int numStateVariables)
44  : G4MagErrorStepper(EqRhs, 6, numStateVariables)
45{
46  // const G4int numberOfVariables= 6;
47  if( numStateVariables < 8 ) 
48  {
49    G4cerr << "ERROR in G4ConstRK4::G4ConstRK4 " 
50           << "   The number of State variables at least 8 " << G4endl;
51    G4cerr << "   Instead it is  numStateVariables= " << numStateVariables << G4endl; 
52    G4Exception("G4ConstRK4::G4ConstRK4()", "InvalidSetup", FatalException,
53        "Valid only for number of state variables of 8 or more. Use another Stepper!");
54  }
55
56  fEq = EqRhs;
57  yMiddle  = new G4double[8];
58  dydxMid  = new G4double[8];
59  yInitial = new G4double[8];
60  yOneStep = new G4double[8];
61
62  dydxm = new G4double[8];
63  dydxt = new G4double[8]; 
64  yt    = new G4double[8]; 
65  Field[0]=0.; Field[1]=0.; Field[2]=0.;
66}
67
68////////////////////////////////////////////////////////////////
69//
70// Destructor
71
72G4ConstRK4::~G4ConstRK4()
73{
74   delete [] yMiddle;
75   delete [] dydxMid;
76   delete [] yInitial;
77   delete [] yOneStep;
78   delete [] dydxm;
79   delete [] dydxt;
80   delete [] yt;
81}
82
83//////////////////////////////////////////////////////////////////////
84//
85// Given values for the variables y[0,..,n-1] and their derivatives
86// dydx[0,...,n-1] known at x, use the classical 4th Runge-Kutta
87// method to advance the solution over an interval h and return the
88// incremented variables as yout[0,...,n-1], which is not a distinct
89// array from y. The user supplies the routine RightHandSide(x,y,dydx),
90// which returns derivatives dydx at x. The source is routine rk4 from
91// NRC p. 712-713 .
92
93void G4ConstRK4::DumbStepper( const G4double  yIn[],
94                              const G4double  dydx[],
95                                    G4double  h,
96                                    G4double  yOut[])
97{
98   G4double  hh = h*0.5 , h6 = h/6.0  ;
99   
100   // 1st Step K1=h*dydx
101   yt[5] = yIn[5] + hh*dydx[5] ;
102   yt[4] = yIn[4] + hh*dydx[4] ;
103   yt[3] = yIn[3] + hh*dydx[3] ;
104   yt[2] = yIn[2] + hh*dydx[2] ;
105   yt[1] = yIn[1] + hh*dydx[1] ;
106   yt[0] = yIn[0] + hh*dydx[0] ;
107   RightHandSideConst(yt,dydxt) ;       
108
109   // 2nd Step K2=h*dydxt
110   yt[5] = yIn[5] + hh*dydxt[5] ;
111   yt[4] = yIn[4] + hh*dydxt[4] ;
112   yt[3] = yIn[3] + hh*dydxt[3] ;
113   yt[2] = yIn[2] + hh*dydxt[2] ;
114   yt[1] = yIn[1] + hh*dydxt[1] ;
115   yt[0] = yIn[0] + hh*dydxt[0] ;
116   RightHandSideConst(yt,dydxm) ;     
117
118   // 3rd Step K3=h*dydxm
119   // now dydxm=(K2+K3)/h
120   yt[5] = yIn[5] + h*dydxm[5] ;
121   dydxm[5] += dydxt[5] ; 
122   yt[4] = yIn[4] + h*dydxm[4] ;
123   dydxm[4] += dydxt[4] ; 
124   yt[3] = yIn[3] + h*dydxm[3] ;
125   dydxm[3] += dydxt[3] ; 
126   yt[2] = yIn[2] + h*dydxm[2] ;
127   dydxm[2] += dydxt[2] ; 
128   yt[1] = yIn[1] + h*dydxm[1] ;
129   dydxm[1] += dydxt[1] ; 
130   yt[0] = yIn[0] + h*dydxm[0] ;
131   dydxm[0] += dydxt[0] ; 
132   RightHandSideConst(yt,dydxt) ;   
133
134   // 4th Step K4=h*dydxt
135   yOut[5] = yIn[5]+h6*(dydx[5]+dydxt[5]+2.0*dydxm[5]);
136   yOut[4] = yIn[4]+h6*(dydx[4]+dydxt[4]+2.0*dydxm[4]);
137   yOut[3] = yIn[3]+h6*(dydx[3]+dydxt[3]+2.0*dydxm[3]);
138   yOut[2] = yIn[2]+h6*(dydx[2]+dydxt[2]+2.0*dydxm[2]);
139   yOut[1] = yIn[1]+h6*(dydx[1]+dydxt[1]+2.0*dydxm[1]);
140   yOut[0] = yIn[0]+h6*(dydx[0]+dydxt[0]+2.0*dydxm[0]);
141   
142}  // end of DumbStepper ....................................................
143
144////////////////////////////////////////////////////////////////
145//
146// Stepper
147
148void
149G4ConstRK4::Stepper( const G4double yInput[],
150                     const G4double dydx[],
151                           G4double hstep,
152                           G4double yOutput[],
153                           G4double yError [] )
154{
155   const G4int nvar = 6;  // number of variables integrated
156   const G4int maxvar= GetNumberOfStateVariables();
157
158   // Correction for Richardson extrapolation
159   G4double  correction = 1. / ( (1 << IntegratorOrder()) -1 );
160
161   G4int i;
162   
163   // Saving yInput because yInput and yOutput can be aliases for same array
164   for (i=0;    i<maxvar; i++) { yInitial[i]= yInput[i]; }
165 
166   // Must copy the part of the state *not* integrated to the output
167   for (i=nvar; i<maxvar; i++) { yOutput[i]=  yInput[i]; }
168
169   // yInitial[7]= yInput[7];  //  The time is typically needed
170   yMiddle[7]  = yInput[7];   // Copy the time from initial value
171   yOneStep[7] = yInput[7];   // As it contributes to final value of yOutput ?
172   // yOutput[7] = yInput[7];  // -> dumb stepper does it too for RK4
173   yError[7] = 0.0;         
174
175   G4double halfStep = hstep * 0.5; 
176
177   // Do two half steps
178   //
179   GetConstField(yInitial,Field);
180   DumbStepper  (yInitial,  dydx,   halfStep, yMiddle);
181   RightHandSideConst(yMiddle, dydxMid);   
182   DumbStepper  (yMiddle, dydxMid, halfStep, yOutput); 
183
184   // Store midpoint, chord calculation
185   //
186   fMidPoint = G4ThreeVector( yMiddle[0],  yMiddle[1],  yMiddle[2]); 
187
188   // Do a full Step
189   //
190   DumbStepper(yInitial, dydx, hstep, yOneStep);
191   for(i=0;i<nvar;i++)
192   {
193      yError [i] = yOutput[i] - yOneStep[i] ;
194      yOutput[i] += yError[i]*correction ;
195        // Provides accuracy increased by 1 order via the
196        // Richardson extrapolation 
197   }
198
199   fInitialPoint = G4ThreeVector( yInitial[0], yInitial[1], yInitial[2]); 
200   fFinalPoint   = G4ThreeVector( yOutput[0],  yOutput[1],  yOutput[2]); 
201
202   return;
203}
204
205////////////////////////////////////////////////////////////////
206//
207// Estimate the maximum distance from the curve to the chord
208//
209// We estimate this using the distance of the midpoint to chord.
210// The method below is good only for angle deviations < 2 pi;
211// this restriction should not be a problem for the Runge Kutta methods,
212// which generally cannot integrate accurately for large angle deviations
213
214G4double G4ConstRK4::DistChord() const 
215{
216  G4double distLine, distChord; 
217
218  if (fInitialPoint != fFinalPoint)
219  {
220     distLine= G4LineSection::Distline( fMidPoint, fInitialPoint, fFinalPoint );
221       // This is a class method that gives distance of Mid
222       // from the Chord between the Initial and Final points
223     distChord = distLine;
224  }
225  else
226  {
227     distChord = (fMidPoint-fInitialPoint).mag();
228  }
229  return distChord;
230}
Note: See TracBrowser for help on using the repository browser.