Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QSStepper.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// G4QSStepper
27//
28// QSS Integrator Stepper
29//
30// Authors: version 1 - Lucio Santi, Rodrigo Castro (Univ. Buenos Aires), 2018-2021
31// version 2 - Mattias Portnoy (Univ. Buenos Aires), 2024
32// --------------------------------------------------------------------
33#ifndef G4QSS_STEPPER_HH
34#define G4QSS_STEPPER_HH
35
36#include "G4FieldTrack.hh"
38#include "G4QSSubstepStruct.hh"
39
40#include <cmath>
42
43/**
44 * @brief G4QSStepper is an integrator of particle's equation of
45 * motion based on the QSS implementation.
46 */
47
49{
50 public:
51
52 /**
53 * Constructor for G4QSStepper.
54 * @param[in] equation Pointer to the provided equation of motion.
55 * @param[in] num_integration_vars The number of integration variables.
56 * @param[in] qssOrder The QSS order (2 or 3 expected; if <= 0 , use value
57 * from Messenger.
58 */
60 G4int num_integration_vars = 6, // always 6 -- ignore
61 G4int qssOrder= -1 );
62
63 /**
64 * Default Destructor. Freeing of memory is done in susbsteps destructor.
65 */
66 ~G4QSStepper() override = default;
67
68 /**
69 * Utility methods.
70 */
71 inline constexpr G4double Cubic_Function(const QSStateVector* states,
72 G4int index, G4double delta_t);
73 inline constexpr G4double Parabolic_Function(const QSStateVector* states,
74 G4int index, G4double delta_t);
75 inline constexpr G4double Linear_Function(const QSStateVector* states,
76 G4int index, G4double delta_t);
77
78 /**
79 * 0 means position type, 1 means velocity type.
80 */
81 inline constexpr int INDEX_TYPE(G4int i);
82
83 /**
84 * Auxiliary methods.
85 */
86 inline void momentum_to_velocity(const G4double* momentum, G4double* out);
87 void set_relativistic_coeff(const G4double* momentum);
89
90 /**
91 * Key methods.
92 */
93 void initialize(const G4double y[]);
94 inline void compare_time_and_update(G4int& index, G4int i);
96 inline void update_field();
98 G4int index, G4double delta_t, G4int order);
100 G4double t, G4double* yOut);
101
102 /**
103 * Moves all the x states of variable index to the current time t.
104 */
105 inline void update_x(G4int index, G4double t);
106
107 /**
108 * Moves all the q states of variable index to the current t.
109 */
110 inline void update_q(G4int index, G4double t);
111
112 /**
113 * Update methods.
114 */
119
120 /*
121 * Updates when does the x,q distance goes beyond the quantum.
122 * Uses polynomial roots-finding formulas.
123 */
124 void update_sync_time(G4int index);
125
126 /**
127 * The stepper for the integration.
128 * The stepsize is fixed, with the step size given by 'h'.
129 * Integrates ODE starting values y[0 to 6]. Outputs yout[].
130 * @param[in] y Starting values array of integration variables.
131 * @param[in] dydx Derivatives array - Not used.
132 * @param[in] h The given step size.
133 * @param[out] yout Integration output.
134 * @param[out] yError The estimated error - Not used.
135 */
136 void Stepper( const G4double y[],
137 const G4double /*dydx*/ [],
138 G4double h,
139 G4double yout[],
140 G4double /* yerr */ [] ) override;
141
142 /**
143 * Returns the QSS order of integration.
144 */
145 inline G4int IntegratorOrder() const override;
146
147 /**
148 * Returns the stepper type-ID, "kQSStepper".
149 */
150 inline G4StepperType StepperType() const override { return kQSStepper; }
151
152 /**
153 * Returns a pointer to the equation of motion.
154 */
156
157 /**
158 * Returns current track state.
159 */
160 inline const field_utils::State& GetYOut() const;
161
162 /**
163 * Track interpolation.
164 * @param[in] tau Step start, x.
165 * @param[in,out] yOut The current track state, y.
166 */
167 void Interpolate(G4double tau, G4double yOut[]);
168
169 /**
170 * Returns the distance from chord line.
171 */
172 inline G4double DistChord() const override;
173
174 /**
175 * Wrapper for the Stepper() function above.
176 */
177 inline void Stepper(const G4double yInput[],
178 const G4double dydx[],
179 G4double hstep,
180 G4double yOutput[],
181 G4double yError[],
182 G4double /*dydxOutput*/ []);
183
184 /**
185 * Sets up interpolation. Does nothing.
186 */
187 inline void SetupInterpolation();
188
189 /*
190 * Obligatory qss driver methods.
191 */
192 inline void reset(const G4FieldTrack* track);
193 inline void SetPrecision(G4double dq_rel, G4double dq_min);
195
196 /*
197 * Sets the mass at rest. Checking/ensuring that it is positive.
198 */
199 inline void setRestMass(G4double restMass);
200
201 private:
202
203 // Constants
204
205 static constexpr int DERIVATIVE_0{0};
206 static constexpr int DERIVATIVE_1{1};
207 static constexpr int DERIVATIVE_2{2};
208 static constexpr int DERIVATIVE_3{3};
209
210 static constexpr int VX{3};
211 static constexpr int VY{4};
212 static constexpr int VZ{5};
213
214 static constexpr int POSITION_IDX{0};
215 static constexpr int VELOCITY_IDX{3};
216
217 static constexpr G4double INFTY{1e+20};
218
219 /** Used to check if field changed from last update field during substeps. */
220 G4bool fField_changed{true};
221 G4bool fTrack_changed{true};
222
223 const G4int qss_order{2};
224
225 Substeps substeps;
226 Substep current_substep;
227 const G4FieldTrack* fCurrent_track{nullptr};
228 QSStateVector dq_vector;
229
230 /** Invariants for this track -- during propagation. */
231 G4double fCharge{-1.0};
232 G4double fCharge_c2;
233 G4double fRestMass{CLHEP::electron_mass_c2};
234 G4double fGamma{1.0};
235 G4double fCoeff; // coeff;
236
237 /** Cached values -- for tiny speed up. */
238 G4double fMassOverC ; // was mass_times_gamma_over_speed_of_light;
239 G4double fInv_mass_over_c;
240
241 /** Used by interpolation driver, need to copy state here when stepper finished. */
242 G4double fYout[12];
243
244 /** QSS parameters separated into velocity and position. */
245 G4double dqrel[2] = {0.0,0.0};
246 G4double dqmin[2] = {0.001,0.001};
247
248 G4double fVelocity{0.0};
249 G4double fFinal_t{0.0};
250};
251
252// ----------------------------------------------------------------------------
253// Inline methods
254// ----------------------------------------------------------------------------
255
256#include "G4QSStepper.icc"
257
258#endif
G4StepperType
G4StepperType defines the available integrator of particle's equation of motion in Geant4.
@ kQSStepper
G4QSStepper.
G4double QSStateVector[6]
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 ...
G4FieldTrack defines a data structure bringing together a magnetic track's state (position,...
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[], G4double[])
G4QSStepper(G4EquationOfMotion *equation, G4int num_integration_vars=6, G4int qssOrder=-1)
void update_x(G4int index, G4double t)
void initialize(const G4double y[])
void update_x_velocity_derivates_using_q(G4int index)
G4int IntegratorOrder() const override
void SetPrecision(G4double dq_rel, G4double dq_min)
void velocity_to_momentum(G4double *y)
void compare_time_and_update(G4int &index, G4int i)
G4int get_next_sync_index()
void update_field()
void Interpolate(G4double tau, G4double yOut[])
void setRestMass(G4double restMass)
void extrapolate_all_states_to_t(Substep *substep, G4double t, G4double *yOut)
void update_x_position_derivates_using_q(G4int index)
G4EquationOfMotion * GetSpecificEquation()
constexpr G4double Cubic_Function(const QSStateVector *states, G4int index, G4double delta_t)
const field_utils::State & GetYOut() const
G4double extrapolate_polynomial(QSStateVector *states, G4int index, G4double delta_t, G4int order)
void update_q(G4int index, G4double t)
void momentum_to_velocity(const G4double *momentum, G4double *out)
G4StepperType StepperType() const override
G4double DistChord() const override
constexpr int INDEX_TYPE(G4int i)
void update_x_derivates_using_q(G4int index)
void set_relativistic_coeff(const G4double *momentum)
void SetupInterpolation()
constexpr G4double Linear_Function(const QSStateVector *states, G4int index, G4double delta_t)
void update_sync_time_one_coefficient(G4int index)
void update_sync_time(G4int index)
G4double GetLastStepLength()
void Stepper(const G4double y[], const G4double[], G4double h, G4double yout[], G4double[]) override
void reset(const G4FieldTrack *track)
~G4QSStepper() override=default
constexpr G4double Parabolic_Function(const QSStateVector *states, G4int index, G4double delta_t)
G4double[G4FieldTrack::ncompSVEC] State
Structs used by G4QSStepper.