Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 // Dormand-Prince RK 6(5) 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: 26 June 2015
32 //
33 // This code is made available subject to the Geant4 license, a copy of
34 // which is available at
35 // http://geant4.org/license
36 // G4DormandPrince745.cc
37 // Geant4
38 //
39 // History
40 // -----------------------------
41 // Created by Somnath on 26 June 2015
42 //
43 //
45 
46 #ifndef DORMAND_PRINCE_RK56_H
47 #define DORMAND_PRINCE_RK56_H
48 
50 
52 
53 
54 public:
55  //constructor
57  G4int numberOfVariables = 6,
58  G4bool primary= true ) ;
59 
60  //destructor
62 
63  //Stepper
64  void Stepper( const G4double y[],
65  const G4double dydx[],
66  G4double h,
67  G4double yout[],
68  G4double yerr[] ) ;
69 
70  G4double DistChord() const;
71  G4int IntegratorOrder() const { return 5; }
72 
75 
76  //For Preparing the Interpolant and calculating the extra stages
77  void SetupInterpolate_low( const G4double yInput[],
78  const G4double dydx[],
79  const G4double Step );
80 
81  //For calculating the output at the tau fraction of Step
82  void Interpolate_low( const G4double yInput[],
83  const G4double dydx[],
84  const G4double Step,
85  G4double yOut[],
86  G4double tau );
87 
89  { SetupInterpolate( fLastInitialVector, fLastDyDx, fLastStepLength); }
90 
91  inline void SetupInterpolate( const G4double yInput[],
92  const G4double dydx[],
93  const G4double Step ){
94  SetupInterpolate_low( yInput, dydx, Step);
95  }
96 
97  //For calculating the output at the tau fraction of Step
98  inline void Interpolate( const G4double yInput[],
99  const G4double dydx[],
100  const G4double Step,
101  G4double yOut[],
102  G4double tau ){
103  Interpolate_low( yInput, dydx, Step, yOut, tau);
104  }
105 
106  void Interpolate( G4double tau, G4double yOut[])
107  { Interpolate( fLastInitialVector, fLastDyDx, fLastStepLength, yOut, tau ); }
108 
109  void SetupInterpolate_high( const G4double yInput[],
110  const G4double dydx[],
111  const G4double Step );
112 
113  //For calculating the output at the tau fraction of Step
114  void Interpolate_high( const G4double yInput[],
115  const G4double dydx[],
116  const G4double Step,
117  G4double yOut[],
118  G4double tau );
119 
120 private:
121  G4double *ak2, *ak3, *ak4, *ak5, *ak6, *ak7, *ak8, *ak9; // for storing intermediate 'k' values in stepper
122  G4double *ak10_low, *ak10, *ak11, * ak12; //For the additional stages of Interpolant
123  G4double *yTemp, *yIn;
124 
125  G4double fLastStepLength;
126  G4double *fLastInitialVector, *fLastFinalVector,
127  *fLastDyDx, *fMidVector, *fMidError;
128  // for DistChord calculations
129 
130  G4DormandPrinceRK56* fAuxStepper;
131 };
132 
133 #endif /* G4DormandPrinceRK56 */
134 
G4double DistChord() const
void Interpolate_high(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)
void Interpolate_low(const G4double yInput[], const G4double dydx[], const G4double Step, G4double yOut[], G4double tau)
void Interpolate(const G4double yInput[], const G4double dydx[], const G4double Step, G4double yOut[], G4double tau)
G4DormandPrinceRK56(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
int G4int
Definition: G4Types.hh:78
void SetupInterpolate_low(const G4double yInput[], const G4double dydx[], const G4double Step)
bool G4bool
Definition: G4Types.hh:79
void Interpolate(G4double tau, G4double yOut[])
Definition: Step.hh:41
G4int IntegratorOrder() const
double G4double
Definition: G4Types.hh:76
G4DormandPrinceRK56 & operator=(const G4DormandPrinceRK56 &)
void SetupInterpolate(const G4double yInput[], const G4double dydx[], const G4double Step)