Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// Author: Somnath Banerjee (CERN, Google Summer of Code 2015), 26.06.2015
33// Supervision: John Apostolakis (CERN)
34// --------------------------------------------------------------------
35#ifndef G4DORMAND_PRINCE_RK56_HH
36#define G4DORMAND_PRINCE_RK56_HH
37
39
40/**
41 * @brief G4DormandPrinceRK56 implements the 6(5) embedded Runge-Kutta
42 * non-FSAL method.
43 */
44
46{
47 public:
48
49 /**
50 * Constructor for G4DormandPrinceRK56.
51 * @param[in] EqRhs Pointer to the provided equation of motion.
52 * @param[in] numberOfVariables The number of integration variables.
53 * @param[in] primary Flag for initialisation of the auxiliary stepper.
54 */
56 G4int numberOfVariables = 6,
57 G4bool primary = true ) ;
58
59 /**
60 * Destructor.
61 */
62 ~G4DormandPrinceRK56() override ;
63
64 /**
65 * Copy constructor and assignment operator not allowed.
66 */
69
70 /**
71 * The stepper for the Runge Kutta integration.
72 * The stepsize is fixed, with the step size given by 'h'.
73 * Integrates ODE starting values y[0 to 6].
74 * Outputs yout[] and its estimated error yerr[].
75 * @param[in] y Starting values array of integration variables.
76 * @param[in] dydx Derivatives array.
77 * @param[in] h The given step size.
78 * @param[out] yout Integration output.
79 * @param[out] yerr The estimated error.
80 */
81 void Stepper( const G4double y[],
82 const G4double dydx[],
83 G4double h,
84 G4double yout[],
85 G4double yerr[] ) override ;
86
87 /**
88 * Returns the distance from chord line.
89 */
90 G4double DistChord() const override;
91
92 /**
93 * Returns the order, 5, of integration.
94 */
95 inline G4int IntegratorOrder() const override { return 5; }
96
97 /**
98 * Returns the stepper type-ID, "kDormandPrinceRK56".
99 */
100 inline G4StepperType StepperType() const override { return kDormandPrinceRK56; }
101
102 /**
103 * Prepares the interpolant and calculates the extra stages.
104 * Fifth order interpolant with one extra function evaluation per step.
105 */
106 void SetupInterpolate_low( const G4double yInput[],
107 const G4double dydx[],
108 const G4double Step );
109
110 /**
111 * Wrappers for SetupInterpolate_low() above.
112 */
113 inline void SetupInterpolate( const G4double yInput[],
114 const G4double dydx[],
115 const G4double Step )
116 {
117 SetupInterpolate_low( yInput, dydx, Step);
118 }
119 inline void SetupInterpolation()
120 {
121 SetupInterpolate( fLastInitialVector, fLastDyDx, fLastStepLength);
122 }
123
124 /**
125 * Calculates the output at the tau fraction of Step.
126 */
127 void Interpolate_low( const G4double yInput[],
128 const G4double dydx[],
129 const G4double Step,
130 G4double yOut[],
131 G4double tau );
132
133 /**
134 * Wrappers for Interpolate_low() above.
135 */
136 inline void Interpolate( const G4double yInput[],
137 const G4double dydx[],
138 const G4double Step,
139 G4double yOut[],
140 G4double tau )
141 {
142 Interpolate_low( yInput, dydx, Step, yOut, tau);
143 }
144 inline void Interpolate( G4double tau, G4double yOut[])
145 {
146 Interpolate( fLastInitialVector, fLastDyDx, fLastStepLength, yOut, tau );
147 }
148
149 /**
150 * Prepares the interpolant and calculates the extra stages.
151 * Sixth order interpolant with 3 additional stages per step.
152 */
153 void SetupInterpolate_high( const G4double yInput[],
154 const G4double dydx[],
155 const G4double Step );
156
157 /**
158 * Calculates the output at the tau fraction of Step, using
159 * the polynomial coefficients and the respective stages.
160 */
161 void Interpolate_high( const G4double yInput[],
162 const G4double dydx[],
163 const G4double Step,
164 G4double yOut[],
165 G4double tau );
166
167 private:
168
169 /** For storing intermediate 'k' values in stepper. */
170 G4double *ak2, *ak3, *ak4, *ak5, *ak6, *ak7, *ak8, *ak9;
171
172 /** For the additional stages of Interpolant. */
173 G4double *ak10_low, *ak10, *ak11, * ak12;
174
175 G4double *yTemp, *yIn;
176
177 G4double fLastStepLength = -1.0;
178
179 /** For DistChord() calculations. */
180 G4double *fLastInitialVector, *fLastFinalVector,
181 *fLastDyDx, *fMidVector, *fMidError;
182
183 G4DormandPrinceRK56* fAuxStepper = nullptr;
184};
185
186#endif
G4StepperType
G4StepperType defines the available integrator of particle's equation of motion in Geant4.
@ kDormandPrinceRK56
G4DormandPrinceRK56.
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4DormandPrinceRK56 implements the 6(5) embedded Runge-Kutta non-FSAL method.
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)
G4int IntegratorOrder() const override
void Interpolate(G4double tau, G4double yOut[])
void Interpolate_low(const G4double yInput[], const G4double dydx[], const G4double Step, G4double yOut[], G4double tau)
void SetupInterpolate(const G4double yInput[], const G4double dydx[], const G4double Step)
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
G4DormandPrinceRK56 & operator=(const G4DormandPrinceRK56 &)=delete
void Interpolate_high(const G4double yInput[], const G4double dydx[], const G4double Step, G4double yOut[], G4double tau)
G4double DistChord() const override
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)
G4StepperType StepperType() const override
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)