Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MagIntegratorDriver.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// G4MagInt_Driver
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 G4MAGINT_DRIVER_HH
37#define G4MAGINT_DRIVER_HH
38
42
43/**
44 * @brief G4MagInt_Driver provides a driver that talks to the Integrator
45 * Stepper and insures that the error is within acceptable bounds.
46 */
47
49 public G4ChordFinderDelegate<G4MagInt_Driver>
50{
51 public:
52
53 /**
54 * Constructor for G4MagInt_Driver.
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 ~G4MagInt_Driver() 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 G4bool DoesReIntegrate() const override { return true; }
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 * Writes out to stream the parameters/state of the driver.
135 */
136 void StreamInfo( std::ostream& os ) const override;
137
138 /**
139 * Attempts one integration step, and returns estimated error 'dyerr'.
140 * It does not ensure accuracy.
141 * @param[in,out] y_posvel The current track in field.
142 * @param[in] dydx dydx array.
143 * @param[in] hstep Proposed step length.
144 * @param[out] dchord_step Estimated sagitta distance.
145 * @param[out] dyerr_pos_sq Estimated error in position.
146 * @param[out] dyerr_mom_rel_sq Estimated error in momentum
147 * (normalised: Delta_Integration(p^2)/(p^2)).
148 * @returns true if integration succeeds.
149 */
150 G4bool QuickAdvance(G4FieldTrack& y_posvel, // In/Out
151 const G4double dydx[],
152 G4double hstep, // In
153 G4double& dchord_step,
154 G4double& dyerr_pos_sq,
155 G4double& dyerr_mom_rel_sq );
156
157 /**
158 * Accessors.
159 */
160 inline G4double GetHmin() const;
161 inline G4double Hmin() const; // Obsolete
162 inline G4double GetSafety() const;
163 inline G4double GetPshrnk() const;
164 inline G4double GetPgrow() const;
165 inline G4double GetErrcon() const;
166 void GetDerivatives(const G4FieldTrack& y_curr, // INput
167 G4double dydx[]) const override; // OUTput
168 void GetDerivatives(const G4FieldTrack& track,
169 G4double dydx[],
170 G4double field[]) const override;
171
172 /**
173 * Getter and setter for the equation of motion.
174 */
176 void SetEquationOfMotion(G4EquationOfMotion* equation) override;
177
178 /**
179 * Sets a new stepper 'pItsStepper' for this driver. Then it calls
180 * ResetParameters() to update its parameters accordingly.
181 */
182 void RenewStepperAndAdjust(G4MagIntegratorStepper* pItsStepper) override;
183
184 /**
185 * Resets the qarameters according to the new provided safety value.
186 * i) sets the exponents (pgrow & pshrnk), using the current order;
187 * ii) sets the safety and calculates "errcon" according to the above values.
188 */
189 inline void ReSetParameters(G4double new_safety = 0.9);
190
191 /**
192 * Modifiers. When setting safety or pgrow, errcon will be set
193 * to a compatible value.
194 */
195 inline void SetSafety(G4double valS);
196 inline void SetPshrnk(G4double valPs);
197 inline void SetPgrow (G4double valPg);
198 inline void SetErrcon(G4double valEc);
200
201 /**
202 * Accessors for the integrator stepper.
203 */
204 const G4MagIntegratorStepper* GetStepper() const override;
206
207 /**
208 * Takes one Step that is as large as possible while satisfying the
209 * accuracy criterion of: yerr < eps * |y_end-y_start|.
210 * @param[in,out] ystart The current track state, y.
211 * @param[in] dydx The derivatives array.
212 * @param[in,out] x Step start, x.
213 * @param[in] htry Step to attempt.
214 * @param[in] eps The relative accuracy.
215 * @param[out] hdid Step achieved.
216 * @param[out] hnext Proposed next step.
217 * @returns true if integration succeeds.
218 */
219 void OneGoodStep(G4double ystart[], // Like old RKF45step()
220 const G4double dydx[],
221 G4double& x,
222 G4double htry,
223 G4double eps,
224 G4double& hdid,
225 G4double& hnext ) ;
226
227 /**
228 * Takes the last step's normalised error and calculates a step size
229 * for the next step. Does it limit the next step's size within a factor
230 * of the current?
231 * -- DOES NOT limit for very bad steps
232 * -- DOES limit for very good (x5).
233 */
234 G4double ComputeNewStepSize(G4double errMaxNorm, // normalised
235 G4double hstepCurrent) override;
236
237 /**
238 * Taking the last step's normalised error, calculates a step size for
239 * the next step. Does not limit the next step's size within a factor of
240 * the current one when *reducing* the size, i.e. for badly failing steps.
241 */
243 G4double hstepCurrent);
244
245 /**
246 * Taking the last step's normalised error, calculates a step size for
247 * the next step. Limits the next step's size within a range around the
248 * current one.
249 */
250 G4double ComputeNewStepSize_WithinLimits(G4double errMaxNorm, // normalised
251 G4double hstepCurrent);
252
253 /**
254 * Modifier and accessor for the maximum number of steps that can be taken
255 * for the integration of a single segment, i.e. a single call to
256 * AccurateAdvance().
257 */
258 inline G4int GetMaxNoSteps() const;
259 inline void SetMaxNoSteps(G4int val);
260
261 /**
262 * More modifiers and accessors.
263 */
264 inline void SetHmin(G4double newval);
265 void SetVerboseLevel(G4int newLevel) override;
266 G4int GetVerboseLevel() const override;
268 void SetSmallestFraction( G4double val );
269
270 protected:
271
272 /**
273 * Loggers, issuing warnings for undesirable situations.
274 */
275 void WarnSmallStepSize(G4double hnext, G4double hstep,
276 G4double h, G4double xDone,
277 G4int noSteps);
278 void WarnTooManyStep(G4double x1start, G4double x2end, G4double xCurrent);
279 void WarnEndPointTooFar(G4double endPointDist,
280 G4double hStepSize ,
281 G4double epsilonRelative,
282 G4int debugFlag);
283
284 /**
285 * Loggers for verbosity printouts.
286 */
287 void PrintStatus(const G4double* StartArr,
288 G4double xstart,
289 const G4double* CurrentArr,
290 G4double xcurrent,
291 G4double requestStep,
292 G4int subStepNo);
293 void PrintStatus(const G4FieldTrack& StartFT,
294 const G4FieldTrack& CurrentFT,
295 G4double requestStep,
296 G4int subStepNo);
297 void PrintStat_Aux(const G4FieldTrack& aFieldTrack,
298 G4double requestStep,
299 G4double actualStep,
300 G4int subStepNo,
301 G4double subStepSize,
302 G4double dotVelocities);
303 /**
304 * Reports on the number of steps, maximum errors etc.
305 */
307
308#ifdef QUICK_ADV_TWO
309 G4bool QuickAdvance( G4double yarrin[], // In
310 const G4double dydx[],
311 G4double hstep,
312 G4double yarrout[], // Out
313 G4double& dchord_step, // Out
314 G4double& dyerr ); // in length
315#endif
316
317 private:
318
319 // ---------------------------------------------------------------
320 // INVARIANTS
321
322 /** Minimum Step allowed in a Step (in absolute units). */
323 G4double fMinimumStep = 0.0;
324
325 /** Smallest fraction of (existing) curve length, in relative units.
326 Below this fraction the current step will be the last. */
327 G4double fSmallestFraction = 1.0e-12; // Expected range 1e-12 to 5e-15
328
329 /** Variables in integration. */
330 const G4int fNoIntegrationVariables = 0;
331
332 /** Minimum number for FieldTrack. */
333 const G4int fMinNoVars = 12;
334
335 /** Full number of variable. */
336 const G4int fNoVars = 0;
337
338 /** Default maximum number of steps is Base divided by the order of Stepper. */
339 G4int fMaxNoSteps;
340 G4int fMaxStepBase = 250; // was 5000
341
342 /** Parameters used to grow and shrink trial stepsize. */
343 G4double safety;
344 G4double pshrnk; // exponent for shrinking
345 G4double pgrow; // exponent for growth
346 G4double errcon;
347
348 G4int fStatisticsVerboseLevel = 0;
349
350 // ---------------------------------------------------------------
351 // DEPENDENT Objects
352
353 G4MagIntegratorStepper* pIntStepper = nullptr;
354
355 // ---------------------------------------------------------------
356 // STATE
357
358 /** Step Statistics. */
359 unsigned long fNoTotalSteps=0, fNoBadSteps=0;
360 unsigned long fNoSmallSteps=0, fNoInitialSmallSteps=0, fNoCalls=0;
361 G4double fDyerr_max=0.0, fDyerr_mx2=0.0;
362 G4double fDyerrPos_smTot=0.0, fDyerrPos_lgTot=0.0, fDyerrVel_lgTot=0.0;
363 G4double fSumH_sm=0.0, fSumH_lg=0.0;
364
365 /** Could be varied during tracking - to help identify issues. */
366 G4int fVerboseLevel = 0; // Verbosity level for printing (debug, ..)
367
368 using ChordFinderDelegate = G4ChordFinderDelegate<G4MagInt_Driver>;
369};
370
371#include "G4MagIntegratorDriver.icc"
372
373#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,...
G4double GetPshrnk() const
void SetPshrnk(G4double valPs)
G4double ComputeNewStepSize_WithinLimits(G4double errMaxNorm, G4double hstepCurrent)
void SetMaxNoSteps(G4int val)
G4MagInt_Driver(const G4MagInt_Driver &)=delete
void SetPgrow(G4double valPg)
void OnStartTracking() override
G4MagInt_Driver(G4double hminimum, G4MagIntegratorStepper *pItsStepper, G4int numberOfComponents=6, G4int statisticsVerbosity=0)
void PrintStatus(const G4double *StartArr, G4double xstart, const G4double *CurrentArr, G4double xcurrent, G4double requestStep, G4int subStepNo)
void SetHmin(G4double newval)
void SetVerboseLevel(G4int newLevel) override
G4double ComputeAndSetErrcon()
void SetEquationOfMotion(G4EquationOfMotion *equation) override
void RenewStepperAndAdjust(G4MagIntegratorStepper *pItsStepper) override
G4double GetErrcon() const
void SetSmallestFraction(G4double val)
const G4MagIntegratorStepper * GetStepper() const override
void SetSafety(G4double valS)
void OnComputeStep(const G4FieldTrack *=nullptr) override
G4MagInt_Driver & operator=(const G4MagInt_Driver &)=delete
G4bool QuickAdvance(G4FieldTrack &y_val, const G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr) override
G4double GetSafety() const
G4double Hmin() const
G4double AdvanceChordLimited(G4FieldTrack &track, G4double stepMax, G4double epsStep, G4double chordDistance) override
G4bool AccurateAdvance(G4FieldTrack &y_current, G4double hstep, G4double eps, G4double hinitial=0.0) override
void WarnSmallStepSize(G4double hnext, G4double hstep, G4double h, G4double xDone, G4int noSteps)
void GetDerivatives(const G4FieldTrack &y_curr, G4double dydx[]) const override
void WarnEndPointTooFar(G4double endPointDist, G4double hStepSize, G4double epsilonRelative, G4int debugFlag)
G4int GetMaxNoSteps() const
void SetErrcon(G4double valEc)
void OneGoodStep(G4double ystart[], const G4double dydx[], G4double &x, G4double htry, G4double eps, G4double &hdid, G4double &hnext)
G4double GetHmin() const
G4double ComputeNewStepSize_WithoutReductionLimit(G4double errMaxNorm, G4double hstepCurrent)
G4int GetVerboseLevel() const override
G4EquationOfMotion * GetEquationOfMotion() override
void PrintStat_Aux(const G4FieldTrack &aFieldTrack, G4double requestStep, G4double actualStep, G4int subStepNo, G4double subStepSize, G4double dotVelocities)
G4double GetSmallestFraction() const
G4bool DoesReIntegrate() const override
G4double GetPgrow() const
G4double ComputeNewStepSize(G4double errMaxNorm, G4double hstepCurrent) override
void WarnTooManyStep(G4double x1start, G4double x2end, G4double xCurrent)
void StreamInfo(std::ostream &os) const override
void ReSetParameters(G4double new_safety=0.9)
G4MagIntegratorStepper is an abstract base class for integrator of particle's equation of motion,...