Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RK547FEq1.hh
Go to the documentation of this file.
1// ********************************************************************
2// * License and Disclaimer *
3// * *
4// * The Geant4 software is copyright of the Copyright Holders of *
5// * the Geant4 Collaboration. It is provided under the terms and *
6// * conditions of the Geant4 Software License, included in the file *
7// * LICENSE and available at http://cern.ch/geant4/license . These *
8// * include a list of copyright holders. *
9// * *
10// * Neither the authors of this software system, nor their employing *
11// * institutes,nor the agencies providing financial support for this *
12// * work make any representation or warranty, express or implied, *
13// * regarding this software system or assume any liability for its *
14// * use. Please see the license in the file LICENSE and URL above *
15// * for the full disclaimer and the limitation of liability. *
16// * *
17// * This code implementation is the result of the scientific and *
18// * technical work of the GEANT4 collaboration. *
19// * By using, copying, modifying or distributing the software (or *
20// * any work based on the software) you agree to acknowledge its *
21// * use in resulting scientific publications, and indicate your *
22// * acceptance of all terms of the Geant4 Software license. *
23// ********************************************************************
24//
25// G4RK547FEq1
26//
27// Class description:
28//
29// An implementation of the 7 stage embedded Runge-Kutta 4,5 pair (RK547FEq1)
30// from the paper:
31// D. J. Higham and G. Hall,
32// "Embedded Runge-Kutta formulae with stable equilibrium states",
33// J. Comput. Appl. Math., vol. 29, no. 1, pp. 25-33, 1990.
34
35// Author: Dmitry Sorokin (CERN, Google Summer of Code 2017), 02.11.2017
36// Supervision: John Apostolakis (CERN)
37// --------------------------------------------------------------------
38#ifndef G4RK547FEQ1_HH
39#define G4RK547FEQ1_HH
40
42#include "G4FieldTrack.hh"
43
44/**
45 * @brief G4RK547FEq1 is an implementation of the 7 stage embedded
46 * Runge-Kutta 4,5 pair.
47 */
48
50{
51 public:
52
53 /**
54 * Constructor for G4RK547FEq1.
55 * @param[in] EqRhs Pointer to the provided equation of motion.
56 * @param[in] integrationVariables The number of integration variables.
57 */
59 G4int integrationVariables = 6);
60
61 /**
62 * Default Destructor.
63 */
64 ~G4RK547FEq1() override = default;
65
66 /**
67 * Copy constructor and assignment operator not allowed.
68 */
69 G4RK547FEq1(const G4RK547FEq1&) = delete;
71
72 /**
73 * The stepper for the Runge Kutta integration.
74 * The stepsize is fixed, with the step size given by 'hstep'.
75 * Integrates ODE starting values yInput[0 to 6].
76 * Outputs yOutput[] and its estimated error yError[].
77 * @param[in] yInput Starting values array of integration variables.
78 * @param[in] dydx Derivatives array.
79 * @param[in] hstep The given step size.
80 * @param[out] yOutput Integration output.
81 * @param[out] yError The estimated error.
82 */
83 void Stepper( const G4double yInput[],
84 const G4double dydx[],
85 G4double hstep,
86 G4double yOutput[],
87 G4double yError[] ) override;
88
89 /**
90 * Same as the Stepper() function above, with dydx also in ouput.
91 * @param[in] yInput Starting values array of integration variables.
92 * @param[in] dydx Derivatives array.
93 * @param[in] hstep The given step size.
94 * @param[out] yOutput Integration output.
95 * @param[out] yError The estimated error.
96 * @param[out] dydxOutput dydx in output.
97 */
98 void Stepper( const G4double yInput[],
99 const G4double dydx[],
100 G4double hstep,
101 G4double yOutput[],
102 G4double yError[],
103 G4double dydxOutput[] );
104
105 /**
106 * Returns the distance from chord line.
107 */
108 G4double DistChord() const override;
109
110 /**
111 * Returns the order, 4, of integration.
112 */
113 inline G4int IntegratorOrder() const override { return 4; }
114
115 /**
116 * Returns the stepper type-ID, "kRK547FEq1".
117 */
118 inline G4StepperType StepperType() const override { return kRK547FEq1; }
119
120 private:
121
122 /**
123 * Utility method used in Stepper() for computing the actual step.
124 */
125 void makeStep( const G4double yInput[],
126 const G4double dydx[],
127 const G4double hstep,
128 G4double yOutput[],
129 G4double* dydxOutput = nullptr,
130 G4double* yError = nullptr) const;
131
132 private:
133
137 fdydxOut[G4FieldTrack::ncompSVEC];
138
139 G4double fhstep = -1.0;
140};
141
142#endif
G4StepperType
G4StepperType defines the available integrator of particle's equation of motion in Geant4.
@ kRK547FEq1
G4RK547FEq1.
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
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)
void Stepper(const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[]) override
G4RK547FEq1(const G4RK547FEq1 &)=delete
G4RK547FEq1 & operator=(const G4RK547FEq1 &)=delete
G4StepperType StepperType() const override
~G4RK547FEq1() override=default
G4RK547FEq1(G4EquationOfMotion *EqRhs, G4int integrationVariables=6)
G4double DistChord() const override
G4int IntegratorOrder() const override