Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CashKarpRKF45.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// G4CashKarpRKF45
27//
28// Class description:
29//
30// The Cash-Karp Runge-Kutta-Fehlberg 4/5 method is an embedded fourth
31// order method (giving fifth-order accuracy) for the solution of an ODE.
32// Two different fourth order estimates are calculated; their difference
33// gives an error estimate. [ref. Numerical Recipes in C, 2nd Edition]
34// It is used to integrate the equations of the motion of a particle
35// in a magnetic field.
36
37// Authors: J.Apostolakis, V.Grichine (CERN), 30.01.1997
38// -------------------------------------------------------------------
39#ifndef G4CASHKARP_RKF45_HH
40#define G4CASHKARP_RKF45_HH
41
43
44/**
45 * @brief G4CashKarpRKF45 implements the Cash-Karp Runge-Kutta-Fehlberg
46 * 4/5 method, an embedded fourth order method (giving fifth-order accuracy)
47 * for the solution of an ODE. Two different fourth order estimates are
48 * calculated; their difference gives an error estimate.
49 * It is used to integrate the equations of the motion of a particle
50 * in a magnetic field.
51 */
52
54{
55 public:
56
57 /**
58 * Constructor for G4CashKarpRKF45.
59 * @param[in] EqRhs Pointer to the provided equation of motion.
60 * @param[in] numberOfVariables The number of integration variables.
61 * @param[in] primary Flag for initialisation of the auxiliary stepper.
62 */
64 G4int numberOfVariables = 6,
65 G4bool primary = true );
66
67 /**
68 * Destructor.
69 */
70 ~G4CashKarpRKF45() override;
71
72 /**
73 * Copy constructor and assignment operator not allowed.
74 */
77
78 /**
79 * The stepper for the Runge Kutta integration.
80 * The stepsize is fixed, with the step size given by 'h'.
81 * Integrates ODE starting values y[0 to 6].
82 * Outputs yout[] and its estimated error yerr[].
83 * @param[in] y Starting values array of integration variables.
84 * @param[in] dydx Derivatives array.
85 * @param[in] h The given step size.
86 * @param[out] yout Integration output.
87 * @param[out] yerr The estimated error.
88 */
89 void Stepper( const G4double y[],
90 const G4double dydx[],
91 G4double h,
92 G4double yout[],
93 G4double yerr[] ) override;
94
95 /**
96 * Returns the distance from chord line.
97 */
98 G4double DistChord() const override;
99
100 /**
101 * Returns the order, 4, of integration.
102 */
103 inline G4int IntegratorOrder() const override { return 4; }
104
105 /**
106 * Returns the stepper type-ID, "kCashKarpRKF45".
107 */
108 inline G4StepperType StepperType() const override { return kCashKarpRKF45; }
109
110 private:
111
112 /** Scratch space. */
113 G4double *ak2, *ak3, *ak4, *ak5, *ak6, *yTemp, *yIn;
114
115 G4double fLastStepLength = 0.0;
116
117 /** For DistChord calculations. */
118 G4double *fLastInitialVector, *fLastFinalVector,
119 *fLastDyDx, *fMidVector, *fMidError;
120
121 G4CashKarpRKF45* fAuxStepper = nullptr;
122};
123
124#endif
G4StepperType
G4StepperType defines the available integrator of particle's equation of motion in Geant4.
@ kCashKarpRKF45
G4CashKarpRKF45.
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4CashKarpRKF45 implements the Cash-Karp Runge-Kutta-Fehlberg 4/5 method, an embedded fourth order me...
~G4CashKarpRKF45() override
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
G4CashKarpRKF45(const G4CashKarpRKF45 &)=delete
G4CashKarpRKF45(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
G4int IntegratorOrder() const override
G4double DistChord() const override
G4CashKarpRKF45 & operator=(const G4CashKarpRKF45 &)=delete
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)