Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DormandPrince745.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// G4DormandPrince745
27//
28// Class desription:
29//
30// An implementation of the 5th order embedded RK method from the paper:
31// J. R. Dormand and P. J. Prince, "A family of embedded Runge-Kutta formulae"
32// Journal of computational and applied Math., vol.6, no.1, pp.19-26, 1980.
33//
34// DormandPrince7 - 5(4) embedded RK method
35
36// Author: Somnath Banerjee (CERN, Google Summer of Code 2015), 25.05.2015
37// Supervision: John Apostolakis (CERN)
38// --------------------------------------------------------------------
39#ifndef G4DORMAND_PRINCE_745_HH
40#define G4DORMAND_PRINCE_745_HH
41
43#include "G4FieldUtils.hh"
44
45/**
46 * @brief G4DormandPrince745 implements the 5th order embedded Runge-Kutta
47 * method, non-FSAL definition of the stepper() method that evaluates one step
48 * in field propagation.
49 */
50
52{
53 public:
54
55 /**
56 * Constructor for G4DormandPrince745.
57 * @param[in] equation Pointer to the provided equation of motion.
58 * @param[in] numberOfVariables The number of integration variables.
59 */
61 G4int numberOfVariables = 6);
62
63 /**
64 * Default Destructor.
65 */
66 ~G4DormandPrince745() override = default;
67
68 /**
69 * Copy constructor and assignment operator not allowed.
70 */
73
74 /**
75 * The stepper for the Runge Kutta integration.
76 * The stepsize is fixed, with the step size given by 'hstep'.
77 * Integrates ODE starting values yInput[0 to 6].
78 * Outputs yOutput[] and its estimated error yError[].
79 * @param[in] yInput Starting values array of integration variables.
80 * @param[in] dydx Derivatives array.
81 * @param[in] hstep The given step size.
82 * @param[out] yOutput Integration output.
83 * @param[out] yError The estimated error.
84 */
85 void Stepper(const G4double yInput[],
86 const G4double dydx[],
87 G4double hstep,
88 G4double yOutput[],
89 G4double yError[]) override;
90
91 /**
92 * Same as the Stepper() function above, with dydx also in ouput.
93 * @param[in] yInput Starting values array of integration variables.
94 * @param[in] dydx Derivatives array.
95 * @param[in] hstep The given step size.
96 * @param[out] yOutput Integration output.
97 * @param[out] yError The estimated error.
98 * @param[out] dydxOutput dysx in output.
99 */
100 void Stepper(const G4double yInput[],
101 const G4double dydx[],
102 G4double hstep,
103 G4double yOutput[],
104 G4double yError[],
105 G4double dydxOutput[]);
106
107
108 /**
109 * Interface method for interpolation setup. Does nothing here.
110 */
111 inline void SetupInterpolation() {}
112
113 /**
114 * Calculates the output at the tau fraction of Step.
115 * Lower (4th) order interpolant given by Dormand and Prince.
116 */
117 void Interpolate4thOrder(G4double yOut[], G4double tau) const;
118
119 /**
120 * Wrapper for Interpolate4thOrder() function above.
121 */
122 inline void Interpolate(G4double tau, G4double yOut[]) const
123 {
124 Interpolate4thOrder(yOut, tau);
125 }
126
127 /**
128 * Sets up the extra stages for the 5th order interpolant.
129 */
131
132 /**
133 * Calculates the interpolated result 'yOut' with the coefficients.
134 * Interpolant of 5th order given by Baker, Dormand, Gilmore and Prince.
135 */
136 void Interpolate5thOrder(G4double yOut[], G4double tau) const;
137
138 /**
139 * Returns the distance from chord line.
140 */
141 G4double DistChord() const override;
142
143 /**
144 * Returns the order, 4, of integration.
145 */
146 inline G4int IntegratorOrder() const override { return 4; }
147
148 /**
149 * Returns the stepper type-ID, "kDormandPrince745".
150 */
151 inline G4StepperType StepperType() const override { return kDormandPrince745; }
152
153 /**
154 * Methods to return the stepper name and description.
155 */
156 const G4String& StepperTypeName() const;
157 const G4String& StepperDescription() const;
158
159 /**
160 * Returns the field state in output.
161 */
162 inline const field_utils::State& GetYOut() const { return fyOut; }
163
164 /**
165 * Returns a pointer to the equation of motion.
166 */
168
169 private:
170
171 field_utils::State ak2, ak3, ak4, ak5, ak6, ak7, ak8, ak9;
172 field_utils::State fyIn, fyOut, fdydxIn;
173
174 G4double fLastStepLength = -1.0;
175};
176
177#endif
G4StepperType
G4StepperType defines the available integrator of particle's equation of motion in Geant4.
@ kDormandPrince745
G4DormandPrince745.
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
void Stepper(const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[]) override
G4int IntegratorOrder() const override
G4StepperType StepperType() const override
G4DormandPrince745(G4EquationOfMotion *equation, G4int numberOfVariables=6)
void Interpolate4thOrder(G4double yOut[], G4double tau) const
G4DormandPrince745(const G4DormandPrince745 &)=delete
G4double DistChord() const override
const G4String & StepperTypeName() const
G4DormandPrince745 & operator=(const G4DormandPrince745 &)=delete
const G4String & StepperDescription() const
~G4DormandPrince745() override=default
const field_utils::State & GetYOut() const
void Interpolate(G4double tau, G4double yOut[]) const
void Interpolate5thOrder(G4double yOut[], G4double tau) const
G4EquationOfMotion * GetSpecificEquation()
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)
G4EquationOfMotion * GetEquationOfMotion()
G4double[G4FieldTrack::ncompSVEC] State