Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ConstRK4.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// G4ConstRK4
27//
28// Class description:
29//
30// G4ConstRK4 performs the integration of one step with error calculation
31// in constant magnetic field. The integration method is the same as in
32// ClassicalRK4. The field value is assumed constant for the step.
33// This field evaluation is called only once per step.
34// G4ConstRK4 can be used only for magnetic fields.
35
36// Authors: J.Apostolakis, T.Nikitina (CERN), 18.09.2008
37// -------------------------------------------------------------------
38#ifndef G4CONSTRK4_HH
39#define G4CONSTRK4_HH
40
41#include "G4MagErrorStepper.hh"
42#include "G4EquationOfMotion.hh"
43#include "G4Mag_EqRhs.hh"
44
45/**
46 * @brief G4ConstRK4 performs the integration of one step with error
47 * calculation in constant magnetic field. The integration method is the
48 * same as in ClassicalRK4. The field value is assumed constant for the step.
49 * This field evaluation is called only once per step.
50 * G4ConstRK4 can be used only for magnetic fields.
51 */
52
54{
55 public:
56
57 /**
58 * Constructor for G4ConstRK4.
59 * @param[in] EqRhs Pointer to the provided equation of motion.
60 * @param[in] numberOfVariables The number of integration variables.
61 */
62 G4ConstRK4(G4Mag_EqRhs* EquationMotion,
63 G4int numberOfStateVariables=8);
64
65 /**
66 * Destructor.
67 */
68 ~G4ConstRK4() override;
69
70 /**
71 * Copy constructor and assignment operator not allowed.
72 */
73 G4ConstRK4(const G4ConstRK4&) = delete;
74 G4ConstRK4& operator=(const G4ConstRK4&) = delete;
75
76 /**
77 * The stepper for the Runge Kutta integration.
78 * The stepsize is fixed, with the step size given by 'h'.
79 * Integrates ODE starting values y[0 to 6].
80 * Outputs yout[] and its estimated error yerr[].
81 * @param[in] y Starting values array of integration variables.
82 * @param[in] dydx Derivatives array.
83 * @param[in] h The given step size.
84 * @param[out] yout Integration output.
85 * @param[out] yerr The estimated error.
86 */
87 void Stepper( const G4double y[],
88 const G4double dydx[],
89 G4double h,
90 G4double yout[],
91 G4double yerr[] ) override;
92
93 /**
94 * Given values for the variables y[0,..,n-1] and their derivatives
95 * dydx[0,...,n-1] known at x, uses the classical 4th Runge-Kutta
96 * method to advance the solution over an interval h and returns the
97 * incremented variables as yout[0,...,n-1]. The user supplies the
98 * function RightHandSide(x,y,dydx), which returns derivatives dydx at x.
99 * The source is routine rk4 from NRC p.712-713.
100 * @param[in] y Starting values array of integration variables.
101 * @param[in] dydx Derivatives array.
102 * @param[in] h The given step size.
103 * @param[out] yout Integration output.
104 */
105 void DumbStepper( const G4double yIn[],
106 const G4double dydx[],
107 G4double h,
108 G4double yOut[] ) override ;
109
110 /**
111 * Returns the distance from chord line.
112 */
113 G4double DistChord() const override;
114
115 /**
116 * Returns the derivatives value, at position and time 'y'.
117 * @param[in] y The position vector plus time (x,y,z,t).
118 * @param[out] dydx The derivatives array.
119 */
120 inline void RightHandSideConst(const G4double y[], G4double dydx[] ) const;
121
122 /**
123 * Returns the field values, at position and time 'y'.
124 * @param[in] y The position vector plus time (x,y,z,t).
125 * @param[out] Field The field value in output.
126 */
127 inline void GetConstField(const G4double y[], G4double Field[]);
128
129 /**
130 * Returns the order, 4, of integration.
131 */
132 inline G4int IntegratorOrder() const override { return 4; }
133
134 /**
135 * Returns the stepper type-ID, "kConstRK4".
136 */
137 inline G4StepperType StepperType() const override { return kConstRK4; }
138
139 private:
140
141 G4ThreeVector fInitialPoint, fMidPoint, fFinalPoint;
142 // Data stored in order to find the chord
143 G4double *dydxm, *dydxt, *yt; // scratch space - not state
144 G4double *yInitial, *yMiddle, *dydxMid, *yOneStep;
145 G4Mag_EqRhs* fEq = nullptr;
146 G4double Field[3];
147};
148
149// Inline methods
150
152 G4double dydx[] ) const
153{
154
155 G4double momentum_mag_square = y[3]*y[3] + y[4]*y[4] + y[5]*y[5];
156 G4double inv_momentum_magnitude = 1.0 / std::sqrt( momentum_mag_square );
157
158 G4double cof = fEq->FCof()*inv_momentum_magnitude;
159
160 dydx[0] = y[3]*inv_momentum_magnitude; // (d/ds)x = Vx/V
161 dydx[1] = y[4]*inv_momentum_magnitude; // (d/ds)y = Vy/V
162 dydx[2] = y[5]*inv_momentum_magnitude; // (d/ds)z = Vz/V
163
164 dydx[3] = cof*(y[4]*Field[2] - y[5]*Field[1]) ; // Ax = a*(Vy*Bz - Vz*By)
165 dydx[4] = cof*(y[5]*Field[0] - y[3]*Field[2]) ; // Ay = a*(Vz*Bx - Vx*Bz)
166 dydx[5] = cof*(y[3]*Field[1] - y[4]*Field[0]) ; // Az = a*(Vx*By - Vy*Bx)
167}
168
170{
171 G4double PositionAndTime[4];
172
173 PositionAndTime[0] = y[0];
174 PositionAndTime[1] = y[1];
175 PositionAndTime[2] = y[2];
176 // Global Time
177 PositionAndTime[3] = y[7];
178 fEq -> GetFieldValue(PositionAndTime, B);
179}
180
181#endif
G4double B(G4double temperature)
G4StepperType
G4StepperType defines the available integrator of particle's equation of motion in Geant4.
@ kConstRK4
G4ConstRK4.
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
G4ConstRK4 & operator=(const G4ConstRK4 &)=delete
G4double DistChord() const override
G4int IntegratorOrder() const override
~G4ConstRK4() override
Definition G4ConstRK4.cc:69
void RightHandSideConst(const G4double y[], G4double dydx[]) const
G4ConstRK4(G4Mag_EqRhs *EquationMotion, G4int numberOfStateVariables=8)
Definition G4ConstRK4.cc:40
void DumbStepper(const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[]) override
Definition G4ConstRK4.cc:90
void GetConstField(const G4double y[], G4double Field[])
G4StepperType StepperType() const override
G4ConstRK4(const G4ConstRK4 &)=delete
G4MagErrorStepper(G4EquationOfMotion *EqRhs, G4int numberOfVariables, G4int numStateVariables=12)
G4Mag_EqRhs is the "standard" equation of motion of a particle in a pure magnetic field.