Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4FSALIntegrationDriver.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: G4FSALIntegrationDriver.hh 109569 2018-05-02 07:08:33Z gcosmo $
28 //
29 //
30 // class G4FSALIntegrationDriver
31 //
32 // Class description:
33 //
34 // Driver class which controls the integration error of a
35 // Runge-Kutta stepper with a FSAL property
36 
37 // History:
38 // - Created. D.Sorokin
39 // --------------------------------------------------------------------
40 
41 #ifndef G4FSALIntegrationDriver_HH
42 #define G4FSALIntegrationDriver_HH
43 
44 #include "G4Types.hh"
45 #include "G4FieldTrack.hh"
46 #include "G4VIntegrationDriver.hh"
47 
48 template <class T>
50 public:
52  G4double hminimum,
53  T* stepper,
54  G4int numberOfComponents = 6,
55  G4int statisticsVerbosity = 1);
56 
57  virtual ~G4FSALIntegrationDriver() override;
58 
61 
62  // Integrates ODE from current s (s=s0) to s=s0+h with accuracy eps.
63  // On output track is replaced by value at end of interval.
64  // The concept is similar to the odeint routine from NRC p.721-722.
65  virtual G4bool AccurateAdvance(
67  G4double hstep,
68  G4double eps, // Requested y_err/hstep
69  G4double hinitial = 0) override; // Suggested 1st interval
70 
71  // QuickAdvance just tries one Step - it does not ensure accuracy.
72  virtual G4bool QuickAdvance(
73  G4FieldTrack& fieldTrack,
74  const G4double dydx[],
75  G4double hstep,
76  G4double& dchord_step,
77  G4double& dyerr) override;
78 
79  virtual void GetDerivatives(
80  const G4FieldTrack &track,
81  G4double dydx[]) const override;
82 
83  // Taking the last step's normalised error, calculate
84  // a step size for the next step.
85  // Do not limit the next step's size within a factor of the
86  // current one.
88  G4double errMaxNorm, // normalised error
89  G4double hstepCurrent) override; // current step size
90 
91  virtual void SetVerboseLevel(G4int newLevel) override;
92  virtual G4int GetVerboseLevel() const override;
93 
94  virtual G4EquationOfMotion* GetEquationOfMotion() override;
95  virtual void SetEquationOfMotion(G4EquationOfMotion* equation) override;
96 
97  virtual const G4MagIntegratorStepper* GetStepper() const override;
98  virtual G4MagIntegratorStepper* GetStepper() override;
99 
100  const T* GetStepperOfPreciseType() const; // Get ptr of precise type
102 
103  // Accessors.
104  G4double GetMinimumStep() const;
105  G4double GetSafety() const;
106  G4double GetPshrnk() const;
107  G4double GetPgrow() const;
108 
109  virtual void RenewStepperAndAdjust(G4MagIntegratorStepper *pItsStepper) override;
110 
111  // Sets a new stepper pItsStepper for this driver. Then it calls
112  // ReSetParameters to reset its parameters accordingly.
113  inline void RenewStepperAndAdjustStrict(T *pItsStepper);
114 
115  // i) sets the exponents (pgrow & pshrnk),
116  // using the current Stepper's order,
117  // ii) sets the safety
118  void ReSetParameters(G4double safety = 0.9);
119 
120 
121  void SetMinimumStep(G4double newval);
122  void SetSafety(G4double valS);
123 
124  // This takes one Step that is as large as possible while
125  // satisfying the accuracy criterion of:
126  // yerr < eps * |y_end-y_start|
127  void OneGoodStep(
128  G4double ystart[],
129  G4double dydx[],
130  G4double& curveLength,
131  G4double htry,
132  G4double eps,
133  G4double& hdid,
134  G4double& hnext);
135 
136  // Modify and Get the Maximum number of Steps that can be
137  // taken for the integration of a single segment -
138  // (ie a single call to AccurateAdvance).
139  G4int GetMaxNoSteps() const;
140  void SetMaxNoSteps( G4int val);
141 
143  void SetSmallestFraction(G4double val);
144 
145 private:
147  G4double GrowStepSize(G4double h, G4double error) const;
148  void UpdateErrorConstraints();
149 
150  void CheckStep(
151  const G4ThreeVector& posIn, const G4ThreeVector& posOut, G4double hdid);
152 
153  // Minimum Step allowed in a Step (in absolute units)
155 
156  // Smallest fraction of (existing) curve length - in relative units
157  // below this fraction the current step will be the last
158  // Expected range 1e-12 to 5e-15;
160 
162 
163  // The (default) maximum number of steps is Base
164  // divided by the order of Stepper
166 
167  // Parameters used to grow and shrink trial stepsize.
169  G4double pshrnk; // exponent for shrinking
170  G4double pgrow; // exponent for growth
171 
172  // muximum error values for shrinking / growing (optimisation).
175 
177 
178  // Step Statistics
180 
181  G4int fVerboseLevel; // Verbosity level for printing (debug, ..)
182  // Could be varied during tracking - to help identify issues
183 
185 };
186 
187 #include "G4FSALIntegrationDriver.icc"
188 
189 #endif
virtual G4double ComputeNewStepSize(G4double errMaxNorm, G4double hstepCurrent) override
G4double GetSafety() const
void CheckStep(const G4ThreeVector &posIn, const G4ThreeVector &posOut, G4double hdid)
G4double GrowStepSize(G4double h, G4double error) const
virtual G4EquationOfMotion * GetEquationOfMotion() override
void ReSetParameters(G4double safety=0.9)
void OneGoodStep(G4double ystart[], G4double dydx[], G4double &curveLength, G4double htry, G4double eps, G4double &hdid, G4double &hnext)
virtual const G4MagIntegratorStepper * GetStepper() const override
G4double GetSmallestFraction() const
void SetSmallestFraction(G4double val)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4int GetMaxNoSteps() const
void SetSafety(G4double valS)
G4double GetPshrnk() const
virtual ~G4FSALIntegrationDriver() override
virtual G4int GetVerboseLevel() const override
const T * GetStepperOfPreciseType() const
void SetMinimumStep(G4double newval)
const G4FSALIntegrationDriver & operator=(const G4FSALIntegrationDriver &)=delete
virtual void SetVerboseLevel(G4int newLevel) override
int G4int
Definition: G4Types.hh:78
virtual void SetEquationOfMotion(G4EquationOfMotion *equation) override
void SetMaxNoSteps(G4int val)
virtual G4bool QuickAdvance(G4FieldTrack &fieldTrack, const G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr) override
static PROLOG_HANDLER error
Definition: xmlrole.cc:112
virtual void GetDerivatives(const G4FieldTrack &track, G4double dydx[]) const override
G4double GetMinimumStep() const
G4double GetPgrow() const
G4double ShrinkStepSize(G4double h, G4double error) const
G4FSALIntegrationDriver(G4double hminimum, T *stepper, G4int numberOfComponents=6, G4int statisticsVerbosity=1)
static const G4double eps
virtual void RenewStepperAndAdjust(G4MagIntegratorStepper *pItsStepper) override
void RenewStepperAndAdjustStrict(T *pItsStepper)
virtual G4bool AccurateAdvance(G4FieldTrack &track, G4double hstep, G4double eps, G4double hinitial=0) override