Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4DormandPrinceRK78.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 // Description of G4DormandPrinceRK78 class
27 //
28 // Implementation of Dormand-Prince 8(7)13M non-FSAL RK method
29 // Author: Somnath Banerjee
30 // Supported by Google as part of Google Summer of Code 2015.
31 // Supervision / code review: John Apostolakis
32 //
33 // Implements a 13 stage embedded explicit Runge-Kutta method,
34 // using a pair of 7th and 8th order formulae
35 //
36 // Paper proposing this RK scheme:
37 // Title: "High order embedded Runge-Kutta formulae",
38 // Authors: P.J. Prince, J.R. Dormand
39 // Journal of Computational and Applied Mathematics, Volume 7, Issue 1, 1981,
40 // Pages 67-75, ISSN 0377-0427,
41 // Reference: DOI: 10.1016/0771-050X(81)90010-3
42 // http://dx.doi.org/10.1016/0771-050X(81)90010-3.
43 // (http://www.sciencedirect.com/science/article/pii/0771050X81900103)
44 //
45 // Created by Somnath on 30/06/15.
46 
47 #ifndef G4Dormand_Prince_RK78_hh
48 #define G4Dormand_Prince_RK78_hh
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  G4double DistChord() const;
67  G4int IntegratorOrder() const {return 7; }
68 
69  private :
72 
73  G4double *ak2, *ak3, *ak4, *ak5, *ak6, *ak7, *ak8,
74  *ak9, *ak10, *ak11, *ak12, *ak13,
75  *yTemp, *yIn;
76 
80  // for DistChord calculations
81 
83 };
84 
85 #endif /* defined(__G4_DormandPrinceRK78_hh) */
G4DormandPrinceRK78 & operator=(const G4DormandPrinceRK78 &)
Float_t y
Definition: compare.C:6
G4int IntegratorOrder() const
G4DormandPrinceRK78(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
G4DormandPrinceRK78 * fAuxStepper
G4double DistChord() const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
int G4int
Definition: G4Types.hh:78