Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ClassicalRK4.cc
Go to the documentation of this file.
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: G4ClassicalRK4.cc 66356 2012-12-18 09:02:32Z gcosmo $
28 //
29 // -------------------------------------------------------------------
30 
31 #include "G4ClassicalRK4.hh"
32 #include "G4ThreeVector.hh"
33 
35 //
36 // Constructor sets the number of variables (default = 6)
37 
39 G4ClassicalRK4(G4EquationOfMotion* EqRhs, G4int numberOfVariables)
40  : G4MagErrorStepper(EqRhs, numberOfVariables)
41 {
42  unsigned int noVariables= std::max(numberOfVariables,8); // For Time .. 7+1
43 
44  dydxm = new G4double[noVariables];
45  dydxt = new G4double[noVariables];
46  yt = new G4double[noVariables];
47 }
48 
50 //
51 // Destructor
52 
54 {
55  delete[] dydxm;
56  delete[] dydxt;
57  delete[] yt;
58 }
59 
61 //
62 // Given values for the variables y[0,..,n-1] and their derivatives
63 // dydx[0,...,n-1] known at x, use the classical 4th Runge-Kutta
64 // method to advance the solution over an interval h and return the
65 // incremented variables as yout[0,...,n-1], which not be a distinct
66 // array from y. The user supplies the routine RightHandSide(x,y,dydx),
67 // which returns derivatives dydx at x. The source is routine rk4 from
68 // NRC p. 712-713 .
69 
70 void
72  const G4double dydx[],
73  G4double h,
74  G4double yOut[])
75 {
76  const G4int nvar = this->GetNumberOfVariables(); // fNumberOfVariables();
77  G4int i;
78  G4double hh = h*0.5 , h6 = h/6.0 ;
79 
80  // Initialise time to t0, needed when it is not updated by the integration.
81  // [ Note: Only for time dependent fields (usually electric)
82  // is it neccessary to integrate the time.]
83  yt[7] = yIn[7];
84  yOut[7] = yIn[7];
85 
86  for(i=0;i<nvar;i++)
87  {
88  yt[i] = yIn[i] + hh*dydx[i] ; // 1st Step K1=h*dydx
89  }
90  RightHandSide(yt,dydxt) ; // 2nd Step K2=h*dydxt
91 
92  for(i=0;i<nvar;i++)
93  {
94  yt[i] = yIn[i] + hh*dydxt[i] ;
95  }
96  RightHandSide(yt,dydxm) ; // 3rd Step K3=h*dydxm
97 
98  for(i=0;i<nvar;i++)
99  {
100  yt[i] = yIn[i] + h*dydxm[i] ;
101  dydxm[i] += dydxt[i] ; // now dydxm=(K2+K3)/h
102  }
103  RightHandSide(yt,dydxt) ; // 4th Step K4=h*dydxt
104 
105  for(i=0;i<nvar;i++) // Final RK4 output
106  {
107  yOut[i] = yIn[i]+h6*(dydx[i]+dydxt[i]+2.0*dydxm[i]); //+K1/6+K4/6+(K2+K3)/3
108  }
109  if ( nvar == 12 ) { NormalisePolarizationVector ( yOut ); }
110 
111 } // end of DumbStepper ....................................................
112 
114 //
115 // StepWithEst
116 
117 void
118 G4ClassicalRK4::StepWithEst( const G4double*,
119  const G4double*,
120  G4double,
121  G4double*,
122  G4double&,
123  G4double&,
124  const G4double*,
125  G4double* )
126 {
127  G4Exception("G4ClassicalRK4::StepWithEst()", "GeomField0001",
128  FatalException, "Method no longer used.");
129 
130 } // end of StepWithEst ......................................................
131 
int G4int
Definition: G4Types.hh:78
G4int GetNumberOfVariables() const
void DumbStepper(const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[])
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void RightHandSide(const double y[], double dydx[])
void NormalisePolarizationVector(G4double vec[12])
G4ClassicalRK4(G4EquationOfMotion *EquationMotion, G4int numberOfVariables=6)
double G4double
Definition: G4Types.hh:76