Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DormandPrinceRK78.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// G4DormandPrinceRK78
27//
28// Class description:
29//
30// Dormand-Prince 8(7)13M non-FSAL RK method, a 13 stage embedded
31// explicit Runge-Kutta method, using a pair of 7th and 8th order formulae.
32//
33// Paper proposing this RK scheme:
34// P.J. Prince, J.R. Dormand, "High order embedded Runge-Kutta formulae",
35// Journal of Computational and Applied Mathematics, Volume 7, Issue 1, 1981,
36// Pages 67-75, ISSN 0377-0427, DOI: 10.1016/0771-050X(81)90010-3
37
38// Author: Somnath Banerjee (CERN, Google Summer of Code 2015), 28.06.2015
39// Supervision: John Apostolakis (CERN)
40// --------------------------------------------------------------------
41#ifndef G4DORMAND_PRINCE_RK78_HH
42#define G4DORMAND_PRINCE_RK78_HH
43
45
46/**
47 * @brief G4DormandPrinceRK78 implements the Dormand-Prince 8(7)13M non-FSAL
48 * Runge-Kutta method, a 13 stage embedded explicit Runge-Kutta method, using
49 * a pair of 7th and 8th order formulae.
50 */
51
53{
54 public:
55
56 /**
57 * Constructor for G4DormandPrince745.
58 * @param[in] EqRhs Pointer to the provided equation of motion.
59 * @param[in] numberOfVariables The number of integration variables.
60 * @param[in] primary Flag for initialisation of the auxiliary stepper.
61 */
63 G4int numberOfVariables = 6,
64 G4bool primary = true);
65
66 /**
67 * Destructor.
68 */
69 ~G4DormandPrinceRK78() override;
70
71 /**
72 * Copy constructor and assignment operator not allowed.
73 */
76
77 /**
78 * The stepper for the Runge Kutta integration.
79 * The stepsize is fixed, with the step size given by 'h'.
80 * Integrates ODE starting values y[0 to 6].
81 * Outputs yout[] and its estimated error yerr[].
82 * @param[in] y Starting values array of integration variables.
83 * @param[in] dydx Derivatives array.
84 * @param[in] h The given step size.
85 * @param[out] yout Integration output.
86 * @param[out] yerr The estimated error.
87 */
88 void Stepper( const G4double y[],
89 const G4double dydx[],
90 G4double h,
91 G4double yout[],
92 G4double yerr[]) override ;
93
94 /**
95 * Returns the distance from chord line.
96 */
97 G4double DistChord() const override;
98
99 /**
100 * Returns the order, 7, of integration.
101 */
102 inline G4int IntegratorOrder() const override { return 7; }
103
104 /**
105 * Returns the stepper type-ID, "kDormandPrinceRK78".
106 */
107 inline G4StepperType StepperType() const override { return kDormandPrinceRK78; }
108
109 private :
110
111 G4double *ak2, *ak3, *ak4, *ak5, *ak6, *ak7, *ak8,
112 *ak9, *ak10, *ak11, *ak12, *ak13,
113 *yTemp, *yIn;
114
115 G4double fLastStepLength = -1.0;
116
117 /** For DistChord() calculations. */
118 G4double *fLastInitialVector, *fLastFinalVector,
119 *fLastDyDx, *fMidVector, *fMidError;
120
121 G4DormandPrinceRK78* fAuxStepper = nullptr;
122};
123
124#endif
G4StepperType
G4StepperType defines the available integrator of particle's equation of motion in Geant4.
@ kDormandPrinceRK78
G4DormandPrinceRK78.
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4DormandPrinceRK78 implements the Dormand-Prince 8(7)13M non-FSAL Runge-Kutta method,...
G4DormandPrinceRK78 & operator=(const G4DormandPrinceRK78 &)=delete
G4DormandPrinceRK78(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
G4StepperType StepperType() const override
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
G4int IntegratorOrder() const override
G4double DistChord() const override
G4DormandPrinceRK78(const G4DormandPrinceRK78 &)=delete
G4EquationOfMotion is the abstract base class for the right hand size of the equation of motion of a ...
G4MagIntegratorStepper(G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, G4bool isFSAL=false)