Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4DormandPrince745.hh
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 // $Id: G4DormandPrince745.hh 97387 2016-06-02 10:03:42Z gcosmo $
27 //
28 // Class desription:
29 // An implementation of the 5th order embedded RK method from the paper
30 // J. R. Dormand and P. J. Prince, “A family of embedded Runge-Kutta formulae,”
31 // Journal of computational and applied …, vol. 6, no. 1, pp. 19–26, 1980.
32 //
33 // DormandPrince7 - 5(4) embedded RK method
34 //
35 // Design & Implementation by Somnath Banerjee
36 // Supervision & code review: John Apostolakis
37 //
38 // Work supported by the Google Summer of Code 2015.
39 //
40 // History
41 // ------------------------------------------
42 // Created : 25 May 2015. - Somnath
43 // Revisions :
44 // * 29 June 2015: Added interpolate() method(s) - Somnath
45 // * May 2016: Cleanup and first comming in G4 - John Apostolakis
46 
47 #ifndef Dormand_Prince_745
48 #define Dormand_Prince_745
49 
51 
53 {
54  public:
56  G4int numberOfVariables = 6,
57  G4bool primary = true);
59 
60  void Stepper( const G4double y[],
61  const G4double dydx[],
62  G4double h,
63  G4double yout[],
64  G4double yerr[] ) ;
65 
66  //For Preparing the Interpolant and calculating the extra stages
67  void SetupInterpolation_low( /* const G4double yInput[],
68  const G4double dydx[],
69  const G4double Step */ );
70 
71  //For calculating the output at the tau fraction of Step
72  void Interpolate_low( /* const G4double yInput[],
73  const G4double dydx[],
74  const G4double Step, */
75  G4double yOut[],
76  G4double tau );
77 
78  inline void SetupInterpolation()
79  /* ( const G4double yInput[],
80  const G4double dydx[],
81  const G4double Step ) */
82  {
83  SetupInterpolation_low( /* yInput, dydx, Step */ );
84  // SetupInterpolation_high( /* yInput, dydx, Step */ );
85  }
86 
87  //For calculating the output at the tau fraction of Step
88  inline void Interpolate(
89  /* const G4double yInput[],
90  const G4double dydx[],
91  const G4double Step, */
92  G4double tau,
93  G4double yOut[]
94  )
95  {
96  Interpolate_low( /* yInput, dydx, Step, */ yOut, tau);
97  // Interpolate_high( /* yInput, dydx, Step, */ yOut, tau);
98  }
99 
100  void SetupInterpolation_high( /* const G4double yInput[],
101  const G4double dydx[],
102  const G4double Step */ );
103 
104  //For calculating the output at the tau fraction of Step
105  void Interpolate_high( /* const G4double yInput[],
106  const G4double dydx[],
107  const G4double Step, */
108  G4double yOut[],
109  G4double tau );
110 
111  G4double DistChord() const;
112  G4double DistChord2() const;
113  G4double DistChord3() const;
114 
115  // Enabling method, with common code between implementations (and steppers)
116  G4double DistLine( G4double yStart[], G4double yMid[], G4double yEnd[] ) const;
117  G4int IntegratorOrder() const {return 4; }
118 
119  //New copy constructor
120  // G4DormandPrince745(const G4DormandPrince745 &);
121 
122 private :
123 
124  G4DormandPrince745& operator=(const G4DormandPrince745&);
125 
126  G4double *ak2, *ak3, *ak4, *ak5, *ak6, *ak7,
127  *ak8, *ak9, //For additional stages in the interpolant
128  *yTemp, *yIn;
129 
130  G4double fLastStepLength;
131  G4double *fLastInitialVector, *fLastFinalVector,
132  *fInitialDyDx, *fMidVector, *fMidError;
133  // for DistChord calculations
134 
135  G4DormandPrince745* fAuxStepper;
136 };
137 #endif /* defined(__Geant4__G4DormandPrince745__) */
G4DormandPrince745(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
void Interpolate(G4double tau, G4double yOut[])
G4double DistLine(G4double yStart[], G4double yMid[], G4double yEnd[]) const
int G4int
Definition: G4Types.hh:78
bool G4bool
Definition: G4Types.hh:79
void Interpolate_high(G4double yOut[], G4double tau)
G4double DistChord3() const
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
G4int IntegratorOrder() const
G4double DistChord2() const
double G4double
Definition: G4Types.hh:76
G4double DistChord() const
void Interpolate_low(G4double yOut[], G4double tau)