Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4BogackiShampine23.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 // Bogacki-Shampine - 4 - 3(2) non-FSAL implementation by Somnath Banerjee
27 // Supervision / code review: John Apostolakis
28 //
29 // Sponsored by Google in Google Summer of Code 2015.
30 //
31 // First version: 20 May 2015
32 //
33 // History
34 // -----------------------------
35 // Created by Somnath Banerjee on 20 May 2015
37 
38 
39 /*
40 
41 This contains the stepper function of the G4BogackiShampine23 class
42 
43 The Bogacki shampine method has the following Butcher's tableau
44 
45 0 |
46 1/2|1/2
47 3/4|0 3/4
48 1 |2/9 1/3 4/9
49 -------------------
50  |2/9 1/3 4/9 0
51  |7/24 1/4 1/3 1/8
52 
53 */
54 
55 #include "G4BogackiShampine23.hh"
56 #include "G4LineSection.hh"
57 
58 using namespace std;
59 
60 //Constructor
62  G4int noIntegrationVariables,
63  G4bool primary)
64  : G4MagIntegratorStepper(EqRhs, noIntegrationVariables),
65  fLastStepLength(0.), fAuxStepper(0)
66 {
67  const G4int numberOfVariables = noIntegrationVariables;
68 
69  ak2 = new G4double[numberOfVariables] ;
70  ak3 = new G4double[numberOfVariables] ;
71  ak4 = new G4double[numberOfVariables] ;
72 
73  pseudoDydx_for_DistChord = new G4double[numberOfVariables];
74 
75  const G4int numStateVars = std::max(noIntegrationVariables,
77 
78  yTemp = new G4double[numberOfVariables] ;
79  yIn = new G4double[numberOfVariables] ;
80 
81  fLastInitialVector = new G4double[numStateVars] ;
82  fLastFinalVector = new G4double[numStateVars] ;
83  fLastDyDx = new G4double[numStateVars];
84 
85  fMidVector = new G4double[numStateVars];
86  fMidError = new G4double[numStateVars];
87  if( primary )
88  {
89  fAuxStepper = new G4BogackiShampine23(EqRhs, numberOfVariables, !primary);
90  }
91 }
92 
93 
94 //Destructor
96 {
97  delete[] ak2;
98  delete[] ak3;
99  delete[] ak4;
100 
101  delete[] yTemp;
102  delete[] yIn;
103 
104  delete[] fLastInitialVector;
105  delete[] fLastFinalVector;
106  delete[] fLastDyDx;
107  delete[] fMidVector;
108  delete[] fMidError;
109 
110  delete fAuxStepper;
111 }
112 
113 //******************************************************************************
114 //
115 // Given values for n = 4 variables yIn[0,...,n-1]
116 // known at x, use the 3rd order Bogacki Shampine method
117 // to advance the solution over an interval Step
118 // and return the incremented variables as yOut[0,...,n-1]. Also
119 // return an estimate of the local truncation error yErr[] using the
120 // embedded 2nd order method. The user supplies routine
121 // RightHandSide(y,dydx), which returns derivatives dydx for y .
122 
123 
124 //******************************************************************************
125 
126 
127 void
129  const G4double DyDx[],
130  G4double Step,
131  G4double yOut[],
132  G4double yErr[])
133 {
134 
135  G4int i;
136 
137  const G4double b21 = 0.5 ,
138  b31 = 0. , b32 = 3.0/4.0 ,
139  b41 = 2.0/9.0, b42 = 1.0/3.0 , b43 = 4.0/9.0;
140 
141 
142  const G4double dc1 = b41 - 7.0/24.0 , dc2 = b42 - 1.0/4.0 ,
143  dc3 = b43 - 1.0/3.0 , dc4 = - 0.125 ;
144 
145 
146 
147  // Initialise time to t0, needed when it is not updated by the integration.
148  // [ Note: Only for time dependent fields (usually electric)
149  // is it neccessary to integrate the time.]
150  yOut[7] = yTemp[7] = yIn[7];
151 
152  const G4int numberOfVariables= this->GetNumberOfVariables(); // The number of variables to be integrated over
153 
154  // Saving yInput because yInput and yOut can be aliases for same array
155 
156  for(i=0;i<numberOfVariables;i++)
157  {
158  yIn[i]=yInput[i];
159  }
160  // RightHandSide(yIn, dydx) ; // 1st Step --Not doing, getting passed
161 
162  for(i=0;i<numberOfVariables;i++)
163  {
164  yTemp[i] = yIn[i] + b21*Step*DyDx[i] ;
165  }
166  RightHandSide(yTemp, ak2) ; // 2nd Step
167 
168  for(i=0;i<numberOfVariables;i++)
169  {
170  yTemp[i] = yIn[i] + Step*(b31*DyDx[i] + b32*ak2[i]) ;
171  }
172  RightHandSide(yTemp, ak3) ; // 3rd Step
173 
174  for(i=0;i<numberOfVariables;i++)
175  {
176  yOut[i] = yIn[i] + Step*(b41*DyDx[i] + b42*ak2[i] + b43*ak3[i]) ;
177  }
178  RightHandSide(yOut, ak4) ; // 4th Step
179 
180  for(i=0;i<numberOfVariables;i++)
181  {
182  // yOut[i] = yIn[i] + Step*(c1*DyDx[i]+ c2*ak2[i] + c3*ak3[i] + c4*ak4[i]);
183 
184  yErr[i] = Step*(dc1*DyDx[i] + dc2*ak2[i] + dc3*ak3[i] +
185  dc4*ak4[i] ) ;
186 
187 
188  // Store Input and Final values, for possible use in calculating chord
189  fLastInitialVector[i] = yIn[i] ;
190  fLastFinalVector[i] = yOut[i];
191  fLastDyDx[i] = DyDx[i];
192  }
193  // NormaliseTangentVector( yOut ); // Not wanted
194 
195  fLastStepLength =Step;
196 
197  return ;
198 }
199 
201 {
202  G4double distLine, distChord;
203  G4ThreeVector initialPoint, finalPoint, midPoint;
204 
205  // Store last initial and final points (they will be overwritten in self-Stepper call!)
206  initialPoint = G4ThreeVector( fLastInitialVector[0],
207  fLastInitialVector[1], fLastInitialVector[2]);
208  finalPoint = G4ThreeVector( fLastFinalVector[0],
209  fLastFinalVector[1], fLastFinalVector[2]);
210 
211  // Do half a step using StepNoErr
212 
213  fAuxStepper->Stepper( fLastInitialVector, fLastDyDx, 0.5 * fLastStepLength,
214  fMidVector, fMidError );
215 
216  midPoint = G4ThreeVector( fMidVector[0], fMidVector[1], fMidVector[2]);
217 
218  // Use stored values of Initial and Endpoint + new Midpoint to evaluate
219  // distance of Chord
220 
221 
222  if (initialPoint != finalPoint)
223  {
224  distLine = G4LineSection::Distline( midPoint, initialPoint, finalPoint );
225  distChord = distLine;
226  }
227  else
228  {
229  distChord = (midPoint-initialPoint).mag();
230  }
231  return distChord;
232 }
233 
234 //------Verified-------
CLHEP::Hep3Vector G4ThreeVector
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
G4BogackiShampine23(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
int G4int
Definition: G4Types.hh:78
G4int GetNumberOfVariables() const
bool G4bool
Definition: G4Types.hh:79
G4double DistChord() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4int GetNumberOfStateVariables() const
void RightHandSide(const double y[], double dydx[])
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
double G4double
Definition: G4Types.hh:76