Geant4-11
G4DormandPrinceRK56.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// G4DormandPrinceRK56
27//
28// Class description:
29//
30// Dormand-Prince RK 6(5) non-FSAL method
31
32// Created: Somnath Banerjee, Google Summer of Code 2015, 26 June 2015
33// Supervision: John Apostolakis, CERN
34// --------------------------------------------------------------------
35#ifndef G4DORMAND_PRINCE_RK56_HH
36#define G4DORMAND_PRINCE_RK56_HH
37
39
41{
42 public:
43
45 G4int numberOfVariables = 6,
46 G4bool primary = true ) ;
47
49
52
53 void Stepper( const G4double y[],
54 const G4double dydx[],
55 G4double h,
56 G4double yout[],
57 G4double yerr[] ) ;
58
59 G4double DistChord() const;
60 G4int IntegratorOrder() const { return 5; }
61
62 void SetupInterpolate_low( const G4double yInput[],
63 const G4double dydx[],
64 const G4double Step );
65 // For preparing the Interpolant and calculating the extra stages
66
67 void Interpolate_low( const G4double yInput[],
68 const G4double dydx[],
69 const G4double Step,
70 G4double yOut[],
71 G4double tau );
72 // For calculating the output at the tau fraction of Step
73
74 inline void SetupInterpolation()
75 {
77 }
78
79 inline void SetupInterpolate( const G4double yInput[],
80 const G4double dydx[],
81 const G4double Step )
82 {
83 SetupInterpolate_low( yInput, dydx, Step);
84 }
85
86 inline void Interpolate( const G4double yInput[],
87 const G4double dydx[],
88 const G4double Step,
89 G4double yOut[],
90 G4double tau )
91 {
92 Interpolate_low( yInput, dydx, Step, yOut, tau);
93 }
94 // For calculating the output at the tau fraction of Step
95
96 inline void Interpolate( G4double tau, G4double yOut[])
97 {
99 }
100
101 void SetupInterpolate_high( const G4double yInput[],
102 const G4double dydx[],
103 const G4double Step );
104
105 void Interpolate_high( const G4double yInput[],
106 const G4double dydx[],
107 const G4double Step,
108 G4double yOut[],
109 G4double tau );
110 // For calculating the output at the tau fraction of Step
111
112 private:
113
114 G4double *ak2, *ak3, *ak4, *ak5, *ak6, *ak7, *ak8, *ak9;
115 // For storing intermediate 'k' values in stepper
117 // For the additional stages of Interpolant
119
123 // For DistChord calculations
124
126};
127
128#endif /* G4DormandPrinceRK56 */
129
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4DormandPrinceRK56(const G4DormandPrinceRK56 &)=delete
void SetupInterpolate_low(const G4double yInput[], const G4double dydx[], const G4double Step)
G4DormandPrinceRK56(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
void Interpolate(G4double tau, G4double yOut[])
void Interpolate_low(const G4double yInput[], const G4double dydx[], const G4double Step, G4double yOut[], G4double tau)
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
void SetupInterpolate(const G4double yInput[], const G4double dydx[], const G4double Step)
G4DormandPrinceRK56 * fAuxStepper
G4DormandPrinceRK56 & operator=(const G4DormandPrinceRK56 &)=delete
void Interpolate_high(const G4double yInput[], const G4double dydx[], const G4double Step, G4double yOut[], G4double tau)
G4double DistChord() const
void Interpolate(const G4double yInput[], const G4double dydx[], const G4double Step, G4double yOut[], G4double tau)
void SetupInterpolate_high(const G4double yInput[], const G4double dydx[], const G4double Step)
G4int IntegratorOrder() const