Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4MagIntegratorDriver.hh
이 파일의 문서화 페이지로 가기
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 //
27 // $Id: G4MagIntegratorDriver.hh 109569 2018-05-02 07:08:33Z gcosmo $
28 //
29 //
30 // class G4MagInt_Driver
31 //
32 // Class description:
33 //
34 // Provides a driver that talks to the Integrator Stepper, and insures that
35 // the error is within acceptable bounds.
36 
37 // History:
38 // - Created. J.Apostolakis.
39 // --------------------------------------------------------------------
40 
41 #ifndef G4MagInt_Driver_Def
42 #define G4MagInt_Driver_Def
43 
44 #include "G4VIntegrationDriver.hh"
46 
48 {
49 public: // with description
50 
51  G4MagInt_Driver(G4double hminimum,
52  G4MagIntegratorStepper* pItsStepper,
53  G4int numberOfComponents = 6,
54  G4int statisticsVerbosity = 1);
55  virtual ~G4MagInt_Driver() override;
56  // Constructor, destructor.
57 
58  G4MagInt_Driver(const G4MagInt_Driver&) = delete;
59  G4MagInt_Driver& operator=(const G4MagInt_Driver&) = delete;
60  // deleted Private copy constructor and assignment operator.
61 
62  virtual G4bool AccurateAdvance(G4FieldTrack& y_current,
63  G4double hstep,
64  G4double eps, // Requested y_err/hstep
65  G4double hinitial = 0.0) override; // Suggested 1st interval
66  // Above drivers for integrator (Runge-Kutta) with stepsize control.
67  // Integrates ODE starting values y_current
68  // from current s (s=s0) to s=s0+h with accuracy eps.
69  // On output ystart is replaced by value at end of interval.
70  // The concept is similar to the odeint routine from NRC p.721-722.
71 
72  virtual G4bool QuickAdvance(G4FieldTrack& y_val, // INOUT
73  const G4double dydx[],
74  G4double hstep, // IN
75  G4double& dchord_step,
76  G4double& dyerr) override;
77  // QuickAdvance just tries one Step - it does not ensure accuracy.
78 
79  G4bool QuickAdvance(G4FieldTrack& y_posvel, // INOUT
80  const G4double dydx[],
81  G4double hstep, // IN
82  G4double& dchord_step,
83  G4double& dyerr_pos_sq,
84  G4double& dyerr_mom_rel_sq ) ;
85  // New QuickAdvance that also just tries one Step
86  // (so also does not ensure accuracy)
87  // but does return the errors in position and
88  // momentum (normalised: Delta_Integration(p^2)/(p^2) )
89 
90  inline G4double GetHmin() const;
91  inline G4double Hmin() const; // Obsolete
92  inline G4double GetSafety() const;
93  inline G4double GetPshrnk() const;
94  inline G4double GetPgrow() const;
95  inline G4double GetErrcon() const;
96  virtual void GetDerivatives(const G4FieldTrack &y_curr, // const, INput
97  G4double dydx[]) const override; // OUTput
98  // Accessors.
99 
100  virtual G4EquationOfMotion* GetEquationOfMotion() override;
101  virtual void SetEquationOfMotion(G4EquationOfMotion* equation) override;
102 
103  virtual void RenewStepperAndAdjust(G4MagIntegratorStepper *pItsStepper) override;
104  // Sets a new stepper pItsStepper for this driver. Then it calls
105  // ReSetParameters to reset its parameters accordingly.
106 
107  inline void ReSetParameters(G4double new_safety= 0.9 );
108  // i) sets the exponents (pgrow & pshrnk),
109  // using the current Stepper's order,
110  // ii) sets the safety
111  // ii) calculates "errcon" according to the above values.
112 
113  inline void SetSafety(G4double valS);
114  inline void SetPshrnk(G4double valPs);
115  inline void SetPgrow (G4double valPg);
116  inline void SetErrcon(G4double valEc);
117  // When setting safety or pgrow, errcon will be set to a
118  // compatible value.
119 
120  inline G4double ComputeAndSetErrcon();
121 
122  virtual const G4MagIntegratorStepper* GetStepper() const override;
123  virtual G4MagIntegratorStepper* GetStepper() override;
124 
125  void OneGoodStep(G4double ystart[], // Like old RKF45step()
126  const G4double dydx[],
127  G4double& x,
128  G4double htry,
129  G4double eps, // memb variables ?
130  G4double& hdid,
131  G4double& hnext ) ;
132  // This takes one Step that is as large as possible while
133  // satisfying the accuracy criterion of:
134  // yerr < eps * |y_end-y_start|
135 
136  virtual G4double ComputeNewStepSize(G4double errMaxNorm, // normalised error
137  G4double hstepCurrent) override; // current step size
138  // Taking the last step's normalised error, calculate
139  // a step size for the next step.
140  // Do not limit the next step's size within a factor of the
141  // current one.
142 
143  G4double ComputeNewStepSize_WithinLimits(G4double errMaxNorm, // normalised error
144  G4double hstepCurrent); // current step size
145  // Taking the last step's normalised error, calculate
146  // a step size for the next step.
147  // Limit the next step's size within a range around the current one.
148 
149  inline G4int GetMaxNoSteps() const;
150  inline void SetMaxNoSteps( G4int val);
151  // Modify and Get the Maximum number of Steps that can be
152  // taken for the integration of a single segment -
153  // (ie a single call to AccurateAdvance).
154 
155 public: // without description
156 
157  inline void SetHmin(G4double newval);
158  virtual void SetVerboseLevel(G4int newLevel) override;
159  virtual G4int GetVerboseLevel() const override;
160 
161  inline G4double GetSmallestFraction() const;
162  void SetSmallestFraction( G4double val );
163 
164 protected: // without description
165  void WarnSmallStepSize(G4double hnext, G4double hstep,
166  G4double h, G4double xDone,
167  G4int noSteps);
168 
169  void WarnTooManyStep(G4double x1start, G4double x2end, G4double xCurrent);
170  void WarnEndPointTooFar(G4double endPointDist,
171  G4double hStepSize ,
172  G4double epsilonRelative,
173  G4int debugFlag);
174  // Issue warnings for undesirable situations
175 
176  void PrintStatus(const G4double* StartArr,
177  G4double xstart,
178  const G4double* CurrentArr,
179  G4double xcurrent,
180  G4double requestStep,
181  G4int subStepNo);
182  void PrintStatus(const G4FieldTrack& StartFT,
183  const G4FieldTrack& CurrentFT,
184  G4double requestStep,
185  G4int subStepNo);
186  void PrintStat_Aux(const G4FieldTrack& aFieldTrack,
187  G4double requestStep,
188  G4double actualStep,
189  G4int subStepNo,
190  G4double subStepSize,
191  G4double dotVelocities);
192  // Verbose output for debugging
193 
194  void PrintStatisticsReport();
195  // Report on the number of steps, maximum errors etc.
196 
197 #ifdef QUICK_ADV_TWO
198  G4bool QuickAdvance( G4double yarrin[], // In
199  const G4double dydx[],
200  G4double hstep,
201  G4double yarrout[], // Out
202  G4double& dchord_step, // Out
203  G4double& dyerr ); // in length
204 #endif
205 
206 private:
207  // ---------------------------------------------------------------
208  // INVARIANTS
209 
211  // Minimum Step allowed in a Step (in absolute units)
212  G4double fSmallestFraction; // Expected range 1e-12 to 5e-15;
213  // Smallest fraction of (existing) curve length - in relative units
214  // below this fraction the current step will be the last
215 
216  const G4int fNoIntegrationVariables; // Number of Variables in integration
217  const G4int fMinNoVars; // Minimum number for FieldTrack
218  const G4int fNoVars; // Full number of variable
219 
222 
224  G4double pshrnk; // exponent for shrinking
225  G4double pgrow; // exponent for growth
227  // Parameters used to grow and shrink trial stepsize.
228 
230 
231  // ---------------------------------------------------------------
232  // DEPENDENT Objects
234 
235  // ---------------------------------------------------------------
236  // STATE
237 
239  unsigned long fNoCalls;
243  // Step Statistics
244 
245  G4int fVerboseLevel; // Verbosity level for printing (debug, ..)
246  // Could be varied during tracking - to help identify issues
247 
248 };
249 
250 #include "G4MagIntegratorDriver.icc"
251 
252 #endif /* G4MagInt_Driver_Def */
Float_t x
Definition: compare.C:6
G4double GetSmallestFraction() const
virtual G4bool AccurateAdvance(G4FieldTrack &y_current, G4double hstep, G4double eps, G4double hinitial=0.0) override
void SetErrcon(G4double valEc)
virtual G4int GetVerboseLevel() const override
virtual G4bool QuickAdvance(G4FieldTrack &y_val, const G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr) override
G4MagInt_Driver & operator=(const G4MagInt_Driver &)=delete
G4double GetPshrnk() const
void SetSmallestFraction(G4double val)
virtual G4double ComputeNewStepSize(G4double errMaxNorm, G4double hstepCurrent) override
G4double ComputeAndSetErrcon()
void SetMaxNoSteps(G4int val)
void WarnSmallStepSize(G4double hnext, G4double hstep, G4double h, G4double xDone, G4int noSteps)
const G4int fNoIntegrationVariables
G4MagInt_Driver(G4double hminimum, G4MagIntegratorStepper *pItsStepper, G4int numberOfComponents=6, G4int statisticsVerbosity=1)
unsigned long fNoSmallSteps
void WarnEndPointTooFar(G4double endPointDist, G4double hStepSize, G4double epsilonRelative, G4int debugFlag)
void SetSafety(G4double valS)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4double GetSafety() const
virtual void GetDerivatives(const G4FieldTrack &y_curr, G4double dydx[]) const override
void SetHmin(G4double newval)
G4MagIntegratorStepper * pIntStepper
G4int GetMaxNoSteps() const
void PrintStatus(const G4double *StartArr, G4double xstart, const G4double *CurrentArr, G4double xcurrent, G4double requestStep, G4int subStepNo)
virtual void SetEquationOfMotion(G4EquationOfMotion *equation) override
virtual G4EquationOfMotion * GetEquationOfMotion() override
G4double ComputeNewStepSize_WithinLimits(G4double errMaxNorm, G4double hstepCurrent)
int G4int
Definition: G4Types.hh:78
void SetPshrnk(G4double valPs)
G4double GetPgrow() const
virtual ~G4MagInt_Driver() override
void ReSetParameters(G4double new_safety=0.9)
G4double GetErrcon() const
void SetPgrow(G4double valPg)
static const G4double eps
unsigned long fNoInitialSmallSteps
unsigned long fNoTotalSteps
void OneGoodStep(G4double ystart[], const G4double dydx[], G4double &x, G4double htry, G4double eps, G4double &hdid, G4double &hnext)
G4double Hmin() const
virtual void RenewStepperAndAdjust(G4MagIntegratorStepper *pItsStepper) override
G4double GetHmin() const
void PrintStat_Aux(const G4FieldTrack &aFieldTrack, G4double requestStep, G4double actualStep, G4int subStepNo, G4double subStepSize, G4double dotVelocities)
void WarnTooManyStep(G4double x1start, G4double x2end, G4double xCurrent)
virtual const G4MagIntegratorStepper * GetStepper() const override
virtual void SetVerboseLevel(G4int newLevel) override