Geant4 10.7.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// V.Grichine, 07.10.1996 - Created
34// W.Wander, 28.01.1998 - Added ability for low order integrators
35// J.Apostolakis, 08.11.2001 - Respect minimum step in AccurateAdvance
36// --------------------------------------------------------------------
37#ifndef G4OLD_MAGINT_DRIVER_HH
38#define G4OLD_MAGINT_DRIVER_HH
39
43
45 public G4ChordFinderDelegate<G4OldMagIntDriver>
46{
47 public: // with description
48
50 G4MagIntegratorStepper* pItsStepper,
51 G4int numberOfComponents = 6,
52 G4int statisticsVerbosity = 0);
53 virtual ~G4OldMagIntDriver() override;
54 // Constructor, destructor.
55
58
60 G4double stepMax,
61 G4double epsStep,
62 G4double chordDistance) override;
63
64 inline virtual void OnStartTracking() override;
65 inline virtual void OnComputeStep() override {};
66 virtual G4bool DoesReIntegrate() const override { return true; }
67
68 virtual G4bool AccurateAdvance(G4FieldTrack& y_current,
69 G4double hstep,
70 G4double eps, // Requested y_err/hstep
71 G4double hinitial = 0.0) override;
72 // Above drivers for integrator (Runge-Kutta) with stepsize control.
73 // Integrates ODE starting values y_current
74 // from current s (s=s0) to s=s0+h with accuracy eps.
75 // On output ystart is replaced by value at end of interval.
76 // The concept is similar to the odeint routine from NRC p.721-722.
77
78 virtual G4bool QuickAdvance(G4FieldTrack& y_val, // INOUT
79 const G4double dydx[],
80 G4double hstep,
81 G4double& dchord_step,
82 G4double& dyerr) override;
83 // QuickAdvance just tries one Step - it does not ensure accuracy.
84
85 G4bool QuickAdvance(G4FieldTrack& y_posvel, // INOUT
86 const G4double dydx[],
87 G4double hstep, // IN
88 G4double& dchord_step,
89 G4double& dyerr_pos_sq,
90 G4double& dyerr_mom_rel_sq );
91 // New QuickAdvance that also just tries one Step
92 // (so also does not ensure accuracy)
93 // but does return the errors in position and
94 // momentum (normalised: Delta_Integration(p^2)/(p^2) )
95
96 inline G4double GetHmin() const;
97 inline G4double Hmin() const; // Obsolete
98 inline G4double GetSafety() const;
99 inline G4double GetPshrnk() const;
100 inline G4double GetPgrow() const;
101 inline G4double GetErrcon() const;
102 virtual void GetDerivatives(const G4FieldTrack& y_curr, // INput
103 G4double dydx[]) const override; // OUTput
104
105 virtual void GetDerivatives(const G4FieldTrack& track,
106 G4double dydx[],
107 G4double field[]) const override;
108 // Accessors
109
110 virtual G4EquationOfMotion* GetEquationOfMotion() override;
111 virtual void SetEquationOfMotion(G4EquationOfMotion* equation) override;
112
113 virtual void RenewStepperAndAdjust(G4MagIntegratorStepper* pItsStepper) override;
114 // Sets a new stepper pItsStepper for this driver. Then it calls
115 // ReSetParameters to reset its parameters accordingly.
116
117 inline void ReSetParameters(G4double new_safety = 0.9);
118 // i) sets the exponents (pgrow & pshrnk),
119 // using the current Stepper's order,
120 // ii) sets the safety
121 // ii) calculates "errcon" according to the above values.
122
123 inline void SetSafety(G4double valS);
124 inline void SetPshrnk(G4double valPs);
125 inline void SetPgrow (G4double valPg);
126 inline void SetErrcon(G4double valEc);
127 // When setting safety or pgrow, errcon will be set to a compatible value.
128
130
131 virtual const G4MagIntegratorStepper* GetStepper() const override;
132 virtual G4MagIntegratorStepper* GetStepper() override;
133
134 void OneGoodStep(G4double ystart[], // Like old RKF45step()
135 const G4double dydx[],
136 G4double& x,
137 G4double htry,
138 G4double eps, // memb variables ?
139 G4double& hdid,
140 G4double& hnext ) ;
141 // This takes one Step that is as large as possible while
142 // satisfying the accuracy criterion of:
143 // yerr < eps * |y_end-y_start|
144
145 virtual G4double ComputeNewStepSize(G4double errMaxNorm, // normalised
146 G4double hstepCurrent) override;
147 // Taking the last step's normalised error, calculate
148 // a step size for the next step.
149 // Do not limit the next step's size within a factor of the
150 // current one.
151
152 void StreamInfo( std::ostream& os ) const override;
153
154 G4double ComputeNewStepSize_WithinLimits(G4double errMaxNorm, // normalised
155 G4double hstepCurrent);
156 // Taking the last step's normalised error, calculate
157 // a step size for the next step.
158 // Limit the next step's size within a range around the current one.
159
160 inline G4int GetMaxNoSteps() const;
161 inline void SetMaxNoSteps(G4int val);
162 // Modify and Get the Maximum number of Steps that can be
163 // taken for the integration of a single segment -
164 // (i.e. a single call to AccurateAdvance).
165
166 public: // without description
167
168 inline void SetHmin(G4double newval);
169 virtual void SetVerboseLevel(G4int newLevel) override;
170 virtual G4int GetVerboseLevel() const override;
171
173 void SetSmallestFraction( G4double val );
174
175 protected: // without description
176
177 void WarnSmallStepSize(G4double hnext, G4double hstep,
178 G4double h, G4double xDone,
179 G4int noSteps);
180
181 void WarnTooManyStep(G4double x1start, G4double x2end, G4double xCurrent);
182 void WarnEndPointTooFar(G4double endPointDist,
183 G4double hStepSize ,
184 G4double epsilonRelative,
185 G4int debugFlag);
186 // Issue warnings for undesirable situations
187
188 void PrintStatus(const G4double* StartArr,
189 G4double xstart,
190 const G4double* CurrentArr,
191 G4double xcurrent,
192 G4double requestStep,
193 G4int subStepNo);
194 void PrintStatus(const G4FieldTrack& StartFT,
195 const G4FieldTrack& CurrentFT,
196 G4double requestStep,
197 G4int subStepNo);
198 void PrintStat_Aux(const G4FieldTrack& aFieldTrack,
199 G4double requestStep,
200 G4double actualStep,
201 G4int subStepNo,
202 G4double subStepSize,
203 G4double dotVelocities);
204 // Verbose output for debugging
205
207 // Report on the number of steps, maximum errors etc.
208
209#ifdef QUICK_ADV_TWO
210 G4bool QuickAdvance( G4double yarrin[], // In
211 const G4double dydx[],
212 G4double hstep,
213 G4double yarrout[], // Out
214 G4double& dchord_step, // Out
215 G4double& dyerr ); // in length
216#endif
217
218 private:
219
220 // ---------------------------------------------------------------
221 // INVARIANTS
222
223 G4double fMinimumStep = 0.0;
224 // Minimum Step allowed in a Step (in absolute units)
225 G4double fSmallestFraction = 1.0e-12; // Expected range 1e-12 to 5e-15
226 // Smallest fraction of (existing) curve length - in relative units
227 // below this fraction the current step will be the last
228
229 const G4int fNoIntegrationVariables = 0; // Variables in integration
230 const G4int fMinNoVars = 12; // Minimum number for FieldTrack
231 const G4int fNoVars = 0; // Full number of variable
232
233 G4int fMaxNoSteps;
234 G4int fMaxStepBase = 250; // was 5000
235 // Default maximum number of steps is Base divided by the order of Stepper
236
237 G4double safety;
238 G4double pshrnk; // exponent for shrinking
239 G4double pgrow; // exponent for growth
240 G4double errcon;
241 // Parameters used to grow and shrink trial stepsize.
242
243 G4int fStatisticsVerboseLevel = 0;
244
245 // ---------------------------------------------------------------
246 // DEPENDENT Objects
247
248 G4MagIntegratorStepper* pIntStepper = nullptr;
249
250 // ---------------------------------------------------------------
251 // STATE
252
253 unsigned long fNoTotalSteps=0, fNoBadSteps=0;
254 unsigned long fNoSmallSteps=0, fNoInitialSmallSteps=0, fNoCalls=0;
255 G4double fDyerr_max=0.0, fDyerr_mx2=0.0;
256 G4double fDyerrPos_smTot=0.0, fDyerrPos_lgTot=0.0, fDyerrVel_lgTot=0.0;
257 G4double fSumH_sm=0.0, fSumH_lg=0.0;
258 // Step Statistics
259
260 G4int fVerboseLevel = 0; // Verbosity level for printing (debug, ..)
261 // Could be varied during tracking - to help identify issues
262
264};
265
266#include "G4OldMagIntDriver.icc"
267
268#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void PrintStat_Aux(const G4FieldTrack &aFieldTrack, G4double requestStep, G4double actualStep, G4int subStepNo, G4double subStepSize, G4double dotVelocities)
virtual G4bool QuickAdvance(G4FieldTrack &y_val, const G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr) override
virtual G4bool DoesReIntegrate() const 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)
virtual void OnStartTracking() override
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
virtual 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)
virtual void SetVerboseLevel(G4int newLevel) override
void SetSafety(G4double valS)
void WarnSmallStepSize(G4double hnext, G4double hstep, G4double h, G4double xDone, G4int noSteps)
void ReSetParameters(G4double new_safety=0.9)
virtual ~G4OldMagIntDriver() override
virtual void OnComputeStep() override
virtual void SetEquationOfMotion(G4EquationOfMotion *equation) override
G4double GetSafety() const
G4double GetErrcon() const
void SetSmallestFraction(G4double val)
void StreamInfo(std::ostream &os) const override
void SetMaxNoSteps(G4int val)
virtual G4double ComputeNewStepSize(G4double errMaxNorm, G4double hstepCurrent) override
void SetErrcon(G4double valEc)
virtual void GetDerivatives(const G4FieldTrack &y_curr, G4double dydx[]) const override
virtual G4EquationOfMotion * GetEquationOfMotion() override
virtual void RenewStepperAndAdjust(G4MagIntegratorStepper *pItsStepper) override
virtual G4int GetVerboseLevel() const override
G4double Hmin() const
G4double GetPshrnk() const
virtual G4double AdvanceChordLimited(G4FieldTrack &track, G4double stepMax, G4double epsStep, G4double chordDistance) override
virtual const G4MagIntegratorStepper * GetStepper() const override
G4double ComputeAndSetErrcon()
G4OldMagIntDriver & operator=(const G4OldMagIntDriver &)=delete
void WarnEndPointTooFar(G4double endPointDist, G4double hStepSize, G4double epsilonRelative, G4int debugFlag)
void WarnTooManyStep(G4double x1start, G4double x2end, G4double xCurrent)