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,v 1.20 2007/05/10 10:10:05 japost Exp $ |
---|
28 | // GEANT4 tag $Name: HEAD $ |
---|
29 | // |
---|
30 | // |
---|
31 | // class G4MagInt_Driver |
---|
32 | // |
---|
33 | // Class description: |
---|
34 | // |
---|
35 | // Provides a driver that talks to the Integrator Stepper, and insures that |
---|
36 | // the error is within acceptable bounds. |
---|
37 | |
---|
38 | // History: |
---|
39 | // - Created. J.Apostolakis. |
---|
40 | // -------------------------------------------------------------------- |
---|
41 | |
---|
42 | #ifndef G4MagInt_Driver_Def |
---|
43 | #define G4MagInt_Driver_Def |
---|
44 | |
---|
45 | #include "G4Types.hh" |
---|
46 | #include "G4FieldTrack.hh" |
---|
47 | #include "G4MagIntegratorStepper.hh" |
---|
48 | |
---|
49 | class G4MagInt_Driver |
---|
50 | { |
---|
51 | public: // with description |
---|
52 | |
---|
53 | G4bool AccurateAdvance(G4FieldTrack& y_current, |
---|
54 | G4double hstep, |
---|
55 | G4double eps, // Requested y_err/hstep |
---|
56 | G4double hinitial=0.0); // Suggested 1st interval |
---|
57 | // Above drivers for integrator (Runge-Kutta) with stepsize control. |
---|
58 | // Integrates ODE starting values y_current |
---|
59 | // from current s (s=s0) to s=s0+h with accuracy eps. |
---|
60 | // On output ystart is replaced by value at end of interval. |
---|
61 | // The concept is similar to the odeint routine from NRC p.721-722. |
---|
62 | |
---|
63 | G4bool QuickAdvance( G4FieldTrack& y_val, // INOUT |
---|
64 | const G4double dydx[], |
---|
65 | G4double hstep, // IN |
---|
66 | G4double& dchord_step, |
---|
67 | G4double& dyerr ) ; |
---|
68 | // QuickAdvance just tries one Step - it does not ensure accuracy. |
---|
69 | |
---|
70 | G4bool QuickAdvance( G4FieldTrack& y_posvel, // INOUT |
---|
71 | const G4double dydx[], |
---|
72 | G4double hstep, // IN |
---|
73 | G4double& dchord_step, |
---|
74 | G4double& dyerr_pos_sq, |
---|
75 | G4double& dyerr_mom_rel_sq ) ; |
---|
76 | // New QuickAdvance that also just tries one Step |
---|
77 | // (so also does not ensure accuracy) |
---|
78 | // but does return the errors in position and |
---|
79 | // momentum (normalised: Delta_Integration(p^2)/(p^2) ) |
---|
80 | |
---|
81 | G4MagInt_Driver( G4double hminimum, |
---|
82 | G4MagIntegratorStepper *pItsStepper, |
---|
83 | G4int numberOfComponents=6, |
---|
84 | G4int statisticsVerbosity=1); |
---|
85 | ~G4MagInt_Driver(); |
---|
86 | // Constructor, destructor. |
---|
87 | |
---|
88 | inline G4double GetHmin() const; |
---|
89 | inline G4double Hmin() const; // Obsolete |
---|
90 | inline G4double GetSafety() const; |
---|
91 | inline G4double GetPshrnk() const; |
---|
92 | inline G4double GetPgrow() const; |
---|
93 | inline G4double GetErrcon() const; |
---|
94 | inline void GetDerivatives( const G4FieldTrack &y_curr, // const, INput |
---|
95 | G4double dydx[] ); // OUTput |
---|
96 | // Accessors. |
---|
97 | |
---|
98 | inline void RenewStepperAndAdjust(G4MagIntegratorStepper *pItsStepper); |
---|
99 | // Sets a new stepper pItsStepper for this driver. Then it calls |
---|
100 | // ReSetParameters to reset its parameters accordingly. |
---|
101 | |
---|
102 | inline void ReSetParameters(G4double new_safety= 0.9 ); |
---|
103 | // i) sets the exponents (pgrow & pshrnk), |
---|
104 | // using the current Stepper's order, |
---|
105 | // ii) sets the safety |
---|
106 | // ii) calculates "errcon" according to the above values. |
---|
107 | |
---|
108 | inline void SetSafety(G4double valS); |
---|
109 | inline void SetPshrnk(G4double valPs); |
---|
110 | inline void SetPgrow (G4double valPg); |
---|
111 | inline void SetErrcon(G4double valEc); |
---|
112 | // When setting safety or pgrow, errcon will be set to a |
---|
113 | // compatible value. |
---|
114 | |
---|
115 | inline G4double ComputeAndSetErrcon(); |
---|
116 | |
---|
117 | inline void SetChargeMomentumMass( G4double particleCharge, |
---|
118 | G4double MomentumXc, |
---|
119 | G4double Mass ); |
---|
120 | // Change them in Equation. particleCharge is in e+ units. |
---|
121 | |
---|
122 | inline const G4MagIntegratorStepper* GetStepper() const; |
---|
123 | |
---|
124 | void OneGoodStep( G4double ystart[], // Like old RKF45step() |
---|
125 | const G4double dydx[], |
---|
126 | G4double& x, |
---|
127 | G4double htry, |
---|
128 | G4double eps, // memb variables ? |
---|
129 | G4double& hdid, |
---|
130 | G4double& hnext ) ; |
---|
131 | // This takes one Step that is as large as possible while |
---|
132 | // satisfying the accuracy criterion of: |
---|
133 | // yerr < eps * |y_end-y_start| |
---|
134 | |
---|
135 | G4double ComputeNewStepSize( G4double errMaxNorm, // normalised error |
---|
136 | G4double hstepCurrent); // current step size |
---|
137 | // Taking the last step's normalised error, calculate |
---|
138 | // a step size for the next step. |
---|
139 | // Do not limit the next step's size within a factor of the |
---|
140 | // current one. |
---|
141 | |
---|
142 | G4double ComputeNewStepSize_WithinLimits( |
---|
143 | 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 | inline void SetVerboseLevel(G4int newLevel); |
---|
159 | inline G4double GetVerboseLevel() const; |
---|
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 | void WarnTooManyStep( G4double x1start, G4double x2end, G4double xCurrent); |
---|
169 | void WarnEndPointTooFar (G4double endPointDist, |
---|
170 | G4double hStepSize , |
---|
171 | G4double epsilonRelative, |
---|
172 | G4int debugFlag); |
---|
173 | // Issue warnings for undesirable situations |
---|
174 | |
---|
175 | void PrintStatus( const G4double* StartArr, |
---|
176 | G4double xstart, |
---|
177 | const G4double* CurrentArr, |
---|
178 | G4double xcurrent, |
---|
179 | G4double requestStep, |
---|
180 | G4int subStepNo ); |
---|
181 | void PrintStatus( const G4FieldTrack& StartFT, |
---|
182 | const G4FieldTrack& CurrentFT, |
---|
183 | G4double requestStep, |
---|
184 | G4int subStepNo ); |
---|
185 | void PrintStat_Aux( const G4FieldTrack& aFieldTrack, |
---|
186 | G4double requestStep, |
---|
187 | G4double actualStep, |
---|
188 | G4int subStepNo, |
---|
189 | G4double subStepSize, |
---|
190 | G4double dotVelocities ); |
---|
191 | // Verbose output for debugging |
---|
192 | |
---|
193 | void PrintStatisticsReport() ; |
---|
194 | // Report on the number of steps, maximum errors etc. |
---|
195 | |
---|
196 | #ifdef QUICK_ADV_TWO |
---|
197 | G4bool QuickAdvance( G4double yarrin[], // In |
---|
198 | const G4double dydx[], |
---|
199 | G4double hstep, |
---|
200 | G4double yarrout[], // Out |
---|
201 | G4double& dchord_step, // Out |
---|
202 | G4double& dyerr ); // in length |
---|
203 | #endif |
---|
204 | |
---|
205 | private: |
---|
206 | |
---|
207 | G4MagInt_Driver(const G4MagInt_Driver&); |
---|
208 | G4MagInt_Driver& operator=(const G4MagInt_Driver&); |
---|
209 | // Private copy constructor and assignment operator. |
---|
210 | |
---|
211 | private: |
---|
212 | |
---|
213 | G4double fMinimumStep; |
---|
214 | // Minimum Step allowed in a Step (in absolute units) |
---|
215 | G4double fSmallestFraction; // Expected range 1e-12 to 5e-15; |
---|
216 | // Smallest fraction of (existing) curve length - in relative units |
---|
217 | // below this fraction the current step will be the last |
---|
218 | |
---|
219 | |
---|
220 | const G4int fNoIntegrationVariables; // Number of Variables in integration |
---|
221 | const G4int fMinNoVars; // Minimum number for FieldTrack |
---|
222 | const G4int fNoVars; // Full number of variable |
---|
223 | |
---|
224 | G4MagIntegratorStepper *pIntStepper; |
---|
225 | G4int fMaxNoSteps; |
---|
226 | static const G4int fMaxStepBase; |
---|
227 | |
---|
228 | G4double safety; |
---|
229 | G4double pshrnk; // exponent for shrinking |
---|
230 | G4double pgrow; // exponent for growth |
---|
231 | G4double errcon; |
---|
232 | // Parameters used to grow and shrink trial stepsize. |
---|
233 | |
---|
234 | static const G4double max_stepping_increase; |
---|
235 | static const G4double max_stepping_decrease; |
---|
236 | // Maximum stepsize increase/decrease factors. |
---|
237 | |
---|
238 | G4int fVerboseLevel; |
---|
239 | // Verbosity level for printing (debug, ..) |
---|
240 | |
---|
241 | G4int fNoTotalSteps, fNoBadSteps, fNoSmallSteps, fNoInitialSmallSteps; |
---|
242 | G4double fDyerr_max, fDyerr_mx2; |
---|
243 | G4double fDyerrPos_smTot, fDyerrVel_smTot, fDyerrPos_lgTot, fDyerrVel_lgTot; |
---|
244 | G4double fSumH_sm, fSumH_lg; |
---|
245 | // Step Statistics |
---|
246 | |
---|
247 | G4int fStatisticsVerboseLevel; |
---|
248 | }; |
---|
249 | |
---|
250 | #include "G4MagIntegratorDriver.icc" |
---|
251 | |
---|
252 | #endif /* G4MagInt_Driver_Def */ |
---|