Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4BogackiShampine45.hh
이 파일의 문서화 페이지로 가기
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 // An implementation of the embedded RK method from the following paper
27 // by P. Bogacki and L. F. Shampine:
28 // “An efficient Runge-Kutta (4,5) pair,”
29 // Comput. Math. with Appl., vol. 32, no. 6, pp. 15–28, Sep. 1996.
30 //
31 // An interpolation method provides the value of an intermediate
32 // point in a step -- if a step was sucessful.
33 //
34 // This version can provide the FSAL property of the method,
35 // which allows the reuse of the last derivative in the next step.
36 // but only by using the additional method GetLastDyDx.
37 // (Alternative interface for simpler use of FSAL is under development.)
38 //
39 // Design & Implementation by Somnath Banerjee
40 // Supervision & code review: John Apostolakis
41 //
42 // Work supported by the Google Summer of Code 2015.
43 //
45 
46 #ifndef Bogacki_Shampine_45
47 #define Bogacki_Shampine_45
48 
50 
52 {
53 public:
55  G4int numberOfVariables = 6,
56  G4bool primary = true);
58 
59  void Stepper( const G4double y[],
60  const G4double dydx[],
61  G4double h,
62  G4double yout[],
63  G4double yerr[] ) ;
64 
65  // This Stepper provides 'dense output'. After a successful
66  // step, it is possible to obtain an estimate of the value
67  // of the function at an intermediate point of the interval.
68  // This requires only two additional evaluations of the
69  // derivative (and thus the field).
70 
71  inline void SetupInterpolation()
72  // const G4double yInput[], const G4double dydx[], const G4double Step )
73  {
74  SetupInterpolationHigh(); // ( yInput, dydx, Step);
75  }
76 
77  //For calculating the output at the tau fraction of Step
78  inline void Interpolate( // const G4double yInput[],
79  // const G4double dydx[],
80  // const G4double Step,
81  G4double tau,
82  G4double yOut[]) // Output value
83  {
84  InterpolateHigh( tau, yOut);
85  // InterpolateHigh( yInput, dydx, Step, yOut, tau);
86  }
87 
89  // Was: ( const G4double yInput[], const G4double, const G4double Step );
90 
91  // For calculating the output at the tau fraction of Step
92  void InterpolateHigh( // const G4double yInput[],
93  // const G4double dydx[],
94  // const G4double Step,
95  G4double tau,
96  G4double yOut[] ) const;
97 
98  G4double DistChord() const;
99  G4int IntegratorOrder() const { return 4; }
100 
101  void GetLastDydx( G4double dyDxLast[] );
102 
103  void PrepareConstants(); // Initialise the values of the bi[][] array
104 
105  private :
106 
109 
110  G4double *ak2, *ak3, *ak4, *ak5, *ak6, *ak7, *ak8,
111  *ak9, *ak10, *ak11, *yTemp, *yIn;
112 
113  G4double *p[6];
114 
118  // for DistChord calculations
119 
120  G4BogackiShampine45* fAuxStepper; // For chord - until interpolation is proven
122 
123  // Class constants
124  static bool fPreparedConstants;
125  static G4double bi[12][7];
126 };
127 
128 #endif /* defined(__Geant4__G4BogackiShampine45__) */
void InterpolateHigh(G4double tau, G4double yOut[]) const
void Interpolate(G4double tau, G4double yOut[])
Float_t y
Definition: compare.C:6
const char * p
Definition: xmltok.h:285
G4double DistChord() const
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
G4BogackiShampine45 * fAuxStepper
G4int IntegratorOrder() const
G4BogackiShampine45 & operator=(const G4BogackiShampine45 &)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void GetLastDydx(G4double dyDxLast[])
static G4double bi[12][7]
G4BogackiShampine45(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
int G4int
Definition: G4Types.hh:78