Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MagIntegratorStepper.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// G4MagIntegratorStepper
27//
28// Class description:
29//
30// Abstract base class for integrator of particle's equation of motion,
31// used in tracking in space dependent magnetic field.
32//
33// A Stepper must integrate over NumberOfVariables elements,
34// and also copy (from input to output) any of NoStateVariables
35// not included in the NumberOfVariables.
36//
37// So it is expected that NoStateVariables >= NumberOfVariables
38
39// Author: John Apostolakis (CERN), 15.01.1997
40// --------------------------------------------------------------------
41#ifndef G4MAGINTEGRATORSTEPPER_HH
42#define G4MAGINTEGRATORSTEPPER_HH
43
44#include "G4Types.hh"
45#include "G4EquationOfMotion.hh"
46#include "G4FieldParameters.hh"
49
51
52/**
53 * @brief G4MagIntegratorStepper is an abstract base class for integrator
54 * of particle's equation of motion, used in tracking in space dependent
55 * magnetic field.
56 */
57
59{
60 public:
61
62 /**
63 * Constructor for G4MagIntegratorStepper.
64 * @param[in] Equation Pointer to the provided equation of motion.
65 * @param[in] numIntegrationVariables The number of integration variables.
66 * @param[in] numStateVariables The number of state variables.
67 * @param[in] isFSAL Flag to indicate if it is an FSAL (First Same As Last)
68 * type driver.
69 */
71 G4int numIntegrationVariables,
72 G4int numStateVariables = 12,
73 G4bool isFSAL = false );
74
75 /**
76 * Default virtual Destructor.
77 */
78 virtual ~G4MagIntegratorStepper() = default;
79
80 /**
81 * Copy constructor and assignment operator not allowed.
82 */
85
86 /**
87 * The stepper for the Runge Kutta integration.
88 * The stepsize is fixed, with the step size given by 'h'.
89 * Integrates ODE starting values y[0 to 6].
90 * Outputs yout[] and its estimated error yerr[].
91 * @param[in] y Starting values array of integration variables.
92 * @param[in] dydx Derivatives array.
93 * @param[in] h The given step size.
94 * @param[out] yout Integration output.
95 * @param[out] yerr The estimated error.
96 */
97 virtual void Stepper( const G4double y[],
98 const G4double dydx[],
99 G4double h,
100 G4double yout[],
101 G4double yerr[] ) = 0;
102
103 /**
104 * Estimates the maximum distance of a chord from the true path
105 * over the segment last integrated.
106 */
107 virtual G4double DistChord() const = 0;
108
109 /**
110 * Simple utility function to (re)normalise 'unit velocity' vector.
111 */
112 inline void NormaliseTangentVector( G4double vec[6] );
113
114 /**
115 * Simple utility function to (re)normalise 'unit spin' vector.
116 */
117 inline void NormalisePolarizationVector( G4double vec[12] );
118
119 /**
120 * Utility method to supply the standard Evaluation of the
121 * Right Hand side of the associated equation.
122 */
123 inline void RightHandSide( const G4double y[], G4double dydx[] ) const;
124
125 /**
126 * Calculates 'dydx' and 'field' at point 'y'.
127 */
128 inline void RightHandSide( const G4double y[],
129 G4double dydx[],
130 G4double field[] ) const;
131
132 /**
133 * Returns the number of variables that the stepper will integrate over.
134 */
136
137 /**
138 * Returns the number of variables of state variables (>= above, integration).
139 */
141
142 /**
143 * Returns the order of the integrator, i.e. its error behaviour is of
144 * the order O(h^order).
145 */
146 virtual G4int IntegratorOrder() const = 0;
147
148 /**
149 * Returns the stepper type ID ('kUserStepper').
150 * This function should be overriden in derived classes.
151 */
152 virtual G4StepperType StepperType() const { return kUserStepper; }
153
154 /**
155 * Replacement method - using new data member.
156 */
158
159 /**
160 * Methods returning the pointer to the associated equation of motion.
161 * As some steppers (e.g. RKG3) require other methods of Eq_Rhs this
162 * function allows for access to them.
163 */
166
167 /**
168 * Setter for the equation of motion.
169 */
170 inline void SetEquationOfMotion(G4EquationOfMotion* newEquation);
171
172 /**
173 * Methods for counting/resetting the number of calls to RHS method(s).
174 */
175 inline unsigned long GetfNoRHSCalls();
176 inline void ResetfNORHSCalls();
177
178 /**
179 * Returns true if the stepper is of FSAL (First Same As Last) type.
180 */
181 inline G4bool IsFSAL() const;
182
183 /**
184 * Returns true if the stepper is of QSS (Quantum State Simulation) type.
185 */
186 inline G4bool isQSS() const;
187 inline void SetIsQSS(G4bool val);
188
189 protected:
190
191 /**
192 * Setters for the integration order and FSAL type.
193 */
194 inline void SetIntegrationOrder(G4int order);
195 inline void SetFSAL(G4bool flag = true);
196
197 private:
198
199 G4EquationOfMotion* fEquation_Rhs = nullptr;
200 const G4int fNoIntegrationVariables = 0; // Variables in integration
201 const G4int fNoStateVariables = 0; // Number required for FieldTrack
202
203 /** Counter for calls to RHS method. */
204 mutable unsigned long fNoRHSCalls = 0UL;
205
206 // Parameters of a RK method -- must be shared by all steppers of a type
207 // -- Invariants for a class
208
209 G4int fIntegrationOrder = -1; // must be set by stepper !!!
210 // All ClassicalRK4 steppers are 4th order
211 G4bool fIsFSAL = false;
212 // Depends on RK method & implementation
213 G4bool fIsQSS = false;
214};
215
216#include "G4MagIntegratorStepper.icc"
217
218#endif
G4StepperType
G4StepperType defines the available integrator of particle's equation of motion in Geant4.
@ kUserStepper
User defined stepper.
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
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 NormalisePolarizationVector(G4double vec[12])
G4bool IsFSAL() const
unsigned long GetfNoRHSCalls()
G4EquationOfMotion * GetEquationOfMotion()
G4int GetNumberOfVariables() const
void RightHandSide(const G4double y[], G4double dydx[]) const
void NormaliseTangentVector(G4double vec[6])
virtual G4double DistChord() const =0
virtual ~G4MagIntegratorStepper()=default
virtual void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])=0
G4MagIntegratorStepper & operator=(const G4MagIntegratorStepper &)=delete
void SetIntegrationOrder(G4int order)
G4bool isQSS() const
void SetEquationOfMotion(G4EquationOfMotion *newEquation)
const G4EquationOfMotion * GetEquationOfMotion() const
void SetIsQSS(G4bool val)
G4MagIntegratorStepper(const G4MagIntegratorStepper &)=delete
void SetFSAL(G4bool flag=true)
virtual G4StepperType StepperType() const
void RightHandSide(const G4double y[], G4double dydx[], G4double field[]) const
G4int GetNumberOfStateVariables() const
virtual G4int IntegratorOrder() const =0