Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4OldMagIntDriver.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// G4OldMagIntDriver
27//
28// Class description:
29//
30// Provides a driver that talks to the Integrator Stepper, and insures that
31// the error is within acceptable bounds.
32
33// Author: Vladimir Grichine (CERN), 07.10.1996 - Created
34// W.Wander (MIT), 28.01.1998 - Added ability for low order integrators
35// --------------------------------------------------------------------
36#ifndef G4OLD_MAGINT_DRIVER_HH
37#define G4OLD_MAGINT_DRIVER_HH
38
42
43/**
44 * @brief G4OldMagIntDriver provides a driver that talks to the Integrator
45 * Stepper and insures that the error is within acceptable bounds.
46 */
47
49 public G4ChordFinderDelegate<G4OldMagIntDriver>
50{
51 public:
52
53 /**
54 * Constructor for G4OldMagIntDriver.
55 * @param[in] hminimum The minumum allowed step.
56 * @param[in] pItsStepper Pointer to the integrator stepper.
57 * @param[in] numberOfComponents The number of integration variables.
58 * @param[in] statisticsVerbosity Flag for verbosity.
59 */
61 G4MagIntegratorStepper* pItsStepper,
62 G4int numberOfComponents = 6,
63 G4int statisticsVerbosity = 0);
64
65 /**
66 * Destructor. Provides statistics if verbosity level is greater than 1.
67 */
68 ~G4OldMagIntDriver() override;
69
70 /**
71 * Copy constructor and assignment operator not allowed.
72 */
75
76 /**
77 * Computes the step to take, based on chord limits.
78 * @param[in,out] track The current track in field.
79 * @param[in] stepMax Proposed maximum step length.
80 * @param[in] epsStep Requested accuracy, y_err/hstep.
81 * @param[in] chordDistance Maximum sagitta distance.
82 * @returns The length of step taken.
83 */
85 G4double stepMax,
86 G4double epsStep,
87 G4double chordDistance) override;
88
89 /**
90 * Dispatch interface method for initialisation/reset of driver.
91 */
92 inline void OnStartTracking() override;
93
94 /**
95 * Dispatch interface method for computing step. Does nothing here.
96 */
97 inline void OnComputeStep(const G4FieldTrack* = nullptr) override;
98
99 /**
100 * The driver implements re-integration, so returns true.
101 */
102 inline G4bool DoesReIntegrate() const override;
103
104 /**
105 * Advances integration accurately by relative accuracy better than 'eps'.
106 * @param[in,out] y_current The current track in field.
107 * @param[in] hstep Proposed step length.
108 * @param[in] eps Requested accuracy, y_err/hstep.
109 * @param[in] hinitial Initial minimum integration step.
110 * @returns true if integration succeeds.
111 */
113 G4double hstep,
114 G4double eps, // Requested y_err/hstep
115 G4double hinitial = 0.0) override;
116
117 /**
118 * Attempts one integration step, and returns estimated error 'dyerr'.
119 * It does not ensure accuracy.
120 * @param[in,out] y_val The current track in field.
121 * @param[in] dydx dydx array.
122 * @param[in] hstep Proposed step length.
123 * @param[out] dchord_step Estimated sagitta distance.
124 * @param[out] dyerr Estimated error.
125 * @returns true if integration succeeds.
126 */
127 G4bool QuickAdvance(G4FieldTrack& y_val, // In/Out
128 const G4double dydx[],
129 G4double hstep,
130 G4double& dchord_step,
131 G4double& dyerr) override;
132
133 /**
134 * Attempts one integration step, and returns estimated error 'dyerr'.
135 * It does not ensure accuracy.
136 * @param[in,out] y_posvel The current track in field.
137 * @param[in] dydx dydx array.
138 * @param[in] hstep Proposed step length.
139 * @param[out] dchord_step Estimated sagitta distance.
140 * @param[out] dyerr_pos_sq Estimated error in position.
141 * @param[out] dyerr_mom_rel_sq Estimated error in momentum
142 * (normalised: Delta_Integration(p^2)/(p^2)).
143 * @returns true if integration succeeds.
144 */
145 G4bool QuickAdvance(G4FieldTrack& y_posvel, // In/Out
146 const G4double dydx[],
147 G4double hstep, // In
148 G4double& dchord_step,
149 G4double& dyerr_pos_sq,
150 G4double& dyerr_mom_rel_sq);
151
152 /**
153 * Accessors.
154 */
155 inline G4double GetHmin() const;
156 inline G4double Hmin() const; // Obsolete
157 inline G4double GetSafety() const;
158 inline G4double GetPshrnk() const;
159 inline G4double GetPgrow() const;
160 inline G4double GetErrcon() const;
161 void GetDerivatives(const G4FieldTrack& y_curr, // INput
162 G4double dydx[]) const override; // OUTput
163 void GetDerivatives(const G4FieldTrack& track,
164 G4double dydx[],
165 G4double field[]) const override;
166 /**
167 * Getter and setter for the equation of motion.
168 */
170 void SetEquationOfMotion(G4EquationOfMotion* equation) override;
171
172 /**
173 * Sets a new stepper 'pItsStepper' for this driver. Then it calls
174 * ResetParameters() to update its parameters accordingly.
175 */
176 void RenewStepperAndAdjust(G4MagIntegratorStepper* pItsStepper) override;
177
178 /**
179 * Resets the qarameters according to the new provided safety value.
180 * i) sets the exponents (pgrow & pshrnk), using the current order;
181 * ii) sets the safety and calculates "errcon" according to the above values.
182 */
183 inline void ReSetParameters(G4double new_safety = 0.9);
184
185 /**
186 * Modifiers. When setting safety or pgrow, errcon will be set
187 * to a compatible value.
188 */
189 inline void SetSafety(G4double valS);
190 inline void SetPshrnk(G4double valPs);
191 inline void SetPgrow (G4double valPg);
192 inline void SetErrcon(G4double valEc);
194
195 /**
196 * Accessors for the integrator stepper.
197 */
198 const G4MagIntegratorStepper* GetStepper() const override;
200
201 /**
202 * Takes one Step that is as large as possible while satisfying the
203 * accuracy criterion of: yerr < eps * |y_end-y_start|.
204 * @param[in,out] ystart The current track state, y.
205 * @param[in] dydx The derivatives array.
206 * @param[in,out] x Step start, x.
207 * @param[in] htry Step to attempt.
208 * @param[in] eps The relative accuracy.
209 * @param[out] hdid Step achieved.
210 * @param[out] hnext Proposed next step.
211 * @returns true if integration succeeds.
212 */
213 void OneGoodStep(G4double ystart[], // Like old RKF45step()
214 const G4double dydx[],
215 G4double& x,
216 G4double htry,
217 G4double eps,
218 G4double& hdid,
219 G4double& hnext) ;
220
221 G4double ComputeNewStepSize(G4double errMaxNorm, // normalised
222 G4double hstepCurrent) override;
223 // Taking the last step's normalised error, calculate
224 // a step size for the next step.
225 // Do not limit the next step's size within a factor of the
226 // current one.
227
228 /**
229 * Writes out to stream the parameters/state of the driver.
230 */
231 void StreamInfo( std::ostream& os ) const override;
232
233 /**
234 * Takes the last step's normalised error and calculates a step size
235 * for the next step. Limits the next step's size within a range around
236 * the current one.
237 */
238 G4double ComputeNewStepSize_WithinLimits(G4double errMaxNorm, // normalised
239 G4double hstepCurrent);
240
241 /**
242 * Modifier and accessor for the maximum number of steps that can be taken
243 * for the integration of a single segment, i.e. a single call to
244 * AccurateAdvance().
245 */
246 inline G4int GetMaxNoSteps() const;
247 inline void SetMaxNoSteps(G4int val);
248
249 /**
250 * More modifiers and accessors.
251 */
252 inline void SetHmin(G4double newval);
253 void SetVerboseLevel(G4int newLevel) override;
254 G4int GetVerboseLevel() const override;
256 void SetSmallestFraction( G4double val );
257
258 protected:
259
260 /**
261 * Loggers, issuing warnings for undesirable situations.
262 */
263 void WarnSmallStepSize(G4double hnext, G4double hstep,
264 G4double h, G4double xDone,
265 G4int noSteps);
266 void WarnTooManyStep(G4double x1start, G4double x2end, G4double xCurrent);
267 void WarnEndPointTooFar(G4double endPointDist,
268 G4double hStepSize ,
269 G4double epsilonRelative,
270 G4int debugFlag);
271
272 /**
273 * Loggers for verbosity printouts.
274 */
275 void PrintStatus(const G4double* StartArr,
276 G4double xstart,
277 const G4double* CurrentArr,
278 G4double xcurrent,
279 G4double requestStep,
280 G4int subStepNo);
281 void PrintStatus(const G4FieldTrack& StartFT,
282 const G4FieldTrack& CurrentFT,
283 G4double requestStep,
284 G4int subStepNo);
285 void PrintStat_Aux(const G4FieldTrack& aFieldTrack,
286 G4double requestStep,
287 G4double actualStep,
288 G4int subStepNo,
289 G4double subStepSize,
290 G4double dotVelocities);
291
292 /**
293 * Reports on the number of steps, maximum errors etc.
294 */
296
297#ifdef QUICK_ADV_TWO
298 G4bool QuickAdvance( G4double yarrin[], // In
299 const G4double dydx[],
300 G4double hstep,
301 G4double yarrout[], // Out
302 G4double& dchord_step, // Out
303 G4double& dyerr ); // in length
304#endif
305
306 private:
307
308 // ---------------------------------------------------------------
309 // INVARIANTS
310
311 /** Minimum Step allowed in a Step (in absolute units). */
312 G4double fMinimumStep = 0.0;
313
314 /** Smallest fraction of (existing) curve length, in relative units.
315 Below this fraction the current step will be the last. */
316 G4double fSmallestFraction = 1.0e-12; // Expected range 1e-12 to 5e-15
317
318 /** Variables in integration. */
319 const G4int fNoIntegrationVariables = 0;
320
321 /** Minimum number for FieldTrack. */
322 const G4int fMinNoVars = 12;
323
324 /** Full number of variable. */
325 const G4int fNoVars = 0;
326
327 /** Default maximum number of steps is Base divided by the order of Stepper. */
328 G4int fMaxNoSteps;
329 G4int fMaxStepBase = 250; // was 5000
330
331 /** Parameters used to grow and shrink trial stepsize. */
332 G4double safety;
333 G4double pshrnk; // exponent for shrinking
334 G4double pgrow; // exponent for growth
335 G4double errcon;
336
337 G4int fStatisticsVerboseLevel = 0;
338
339 // ---------------------------------------------------------------
340 // DEPENDENT Objects
341
342 G4MagIntegratorStepper* pIntStepper = nullptr;
343
344 // ---------------------------------------------------------------
345 // STATE
346
347 /** Step Statistics. */
348 unsigned long fNoTotalSteps=0, fNoBadSteps=0;
349 unsigned long fNoSmallSteps=0, fNoInitialSmallSteps=0, fNoCalls=0;
350 G4double fDyerr_max=0.0, fDyerr_mx2=0.0;
351 G4double fDyerrPos_smTot=0.0, fDyerrPos_lgTot=0.0, fDyerrVel_lgTot=0.0;
352 G4double fSumH_sm=0.0, fSumH_lg=0.0;
353
354 /** Could be varied during tracking - to help identify issues. */
355 G4int fVerboseLevel = 0; // Verbosity level for printing (debug, ..)
356
357 using ChordFinderDelegate = G4ChordFinderDelegate<G4OldMagIntDriver>;
358};
359
360#include "G4OldMagIntDriver.icc"
361
362#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4ChordFinderDelegate is a templated class for a common algorithm of finding step size with distance ...
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 is an abstract base class for integrator of particle's equation of motion,...
void PrintStat_Aux(const G4FieldTrack &aFieldTrack, G4double requestStep, G4double actualStep, G4int subStepNo, G4double subStepSize, G4double dotVelocities)
G4bool QuickAdvance(G4FieldTrack &y_val, const G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr) override
void OnComputeStep(const G4FieldTrack *=nullptr) override
void SetPshrnk(G4double valPs)
G4double ComputeNewStepSize_WithinLimits(G4double errMaxNorm, G4double hstepCurrent)
G4int GetMaxNoSteps() const
G4double GetHmin() const
void SetHmin(G4double newval)
void SetPgrow(G4double valPg)
G4OldMagIntDriver(const G4OldMagIntDriver &)=delete
void PrintStatus(const G4double *StartArr, G4double xstart, const G4double *CurrentArr, G4double xcurrent, G4double requestStep, G4int subStepNo)
G4double GetSmallestFraction() const
G4double GetPgrow() const
G4bool AccurateAdvance(G4FieldTrack &y_current, G4double hstep, G4double eps, G4double hinitial=0.0) override
void OneGoodStep(G4double ystart[], const G4double dydx[], G4double &x, G4double htry, G4double eps, G4double &hdid, G4double &hnext)
void SetSafety(G4double valS)
void WarnSmallStepSize(G4double hnext, G4double hstep, G4double h, G4double xDone, G4int noSteps)
void ReSetParameters(G4double new_safety=0.9)
~G4OldMagIntDriver() override
void SetEquationOfMotion(G4EquationOfMotion *equation) override
G4double GetSafety() const
G4double GetErrcon() const
void SetSmallestFraction(G4double val)
void StreamInfo(std::ostream &os) const override
void GetDerivatives(const G4FieldTrack &y_curr, G4double dydx[]) const override
void SetMaxNoSteps(G4int val)
G4double ComputeNewStepSize(G4double errMaxNorm, G4double hstepCurrent) override
void SetErrcon(G4double valEc)
void OnStartTracking() override
G4EquationOfMotion * GetEquationOfMotion() override
void RenewStepperAndAdjust(G4MagIntegratorStepper *pItsStepper) override
G4bool DoesReIntegrate() const override
G4double AdvanceChordLimited(G4FieldTrack &track, G4double stepMax, G4double epsStep, G4double chordDistance) override
G4double Hmin() const
G4double GetPshrnk() const
G4int GetVerboseLevel() const override
const G4MagIntegratorStepper * GetStepper() const override
G4OldMagIntDriver(G4double hminimum, G4MagIntegratorStepper *pItsStepper, G4int numberOfComponents=6, G4int statisticsVerbosity=0)
G4double ComputeAndSetErrcon()
G4OldMagIntDriver & operator=(const G4OldMagIntDriver &)=delete
void SetVerboseLevel(G4int newLevel) override
void WarnEndPointTooFar(G4double endPointDist, G4double hStepSize, G4double epsilonRelative, G4int debugFlag)
void WarnTooManyStep(G4double x1start, G4double x2end, G4double xCurrent)