Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4FSALIntegrationDriver.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 //
27 // $Id: G4FSALIntegrationDriver.hh 97387 2016-06-02 10:03:42Z gcosmo $
28 //
29 //
30 // class G4FSALIntegrationDriver
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 G4FSALIntegrationDriver_Def
42 #define G4FSALIntegrationDriver_Def
43 
44 #include "G4Types.hh"
45 #include "G4FieldTrack.hh"
47 
49 {
50  public: // with description
51 
53  G4double hstep,
54  G4double eps, // Requested y_err/hstep
55  G4double hinitial=0.0); // Suggested 1st interval
56  // Above drivers for integrator (Runge-Kutta) with stepsize control.
57  // Integrates ODE starting values y_current
58  // from current s (s=s0) to s=s0+h with accuracy eps.
59  // On output ystart is replaced by value at end of interval.
60  // The concept is similar to the odeint routine from NRC p.721-722.
61 
62  G4bool QuickAdvance( G4FieldTrack& y_val, // INOUT
63  G4double dydx[],
64  G4double hstep, // IN
65  G4double& dchord_step,
66  G4double& dyerr ) ;
67  // QuickAdvance just tries one Step - it does not ensure accuracy.
68 
69  G4bool QuickAdvance( G4FieldTrack& y_posvel, // INOUT
70  G4double dydx[],
71  G4double hstep, // IN
72  G4double& dchord_step,
73  G4double& dyerr_pos_sq,
74  G4double& dyerr_mom_rel_sq ) ;
75  // New QuickAdvance that also just tries one Step
76  // (so also does not ensure accuracy)
77  // but does return the errors in position and
78  // momentum (normalised: Delta_Integration(p^2)/(p^2) )
79 
81  G4VFSALIntegrationStepper *pItsStepper,
82  G4int numberOfComponents=6,
83  G4int statisticsVerbosity=1);
85  // Constructor, destructor.
86 
87  inline G4double GetHmin() const;
88  inline G4double Hmin() const; // Obsolete
89  inline G4double GetSafety() const;
90  inline G4double GetPshrnk() const;
91  inline G4double GetPgrow() const;
92  inline G4double GetErrcon() const;
93  inline G4int GetNoTotalSteps() const; //Only for debug purposes
94  inline void GetDerivatives( const G4FieldTrack &y_curr, // const, INput
95  G4double dydx[] ); // OUTput
96  // Accessors.
97 
98  inline void RenewStepperAndAdjust(G4VFSALIntegrationStepper *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 const G4VFSALIntegrationStepper* GetStepper() const;
119 
120  void OneGoodStep( G4double ystart[], // Like old RKF45step()
121  G4double dydx[],
122  G4double& x,
123  G4double htry,
124  G4double eps, // memb variables ?
125  G4double& hdid,
126  G4double& hnext ) ;
127  // This takes one Step that is as large as possible while
128  // satisfying the accuracy criterion of:
129  // yerr < eps * |y_end-y_start|
130 
131  G4double ComputeNewStepSize( G4double errMaxNorm, // normalised error
132  G4double hstepCurrent); // current step size
133  // Taking the last step's normalised error, calculate
134  // a step size for the next step.
135  // Do not limit the next step's size within a factor of the
136  // current one.
137 
139  G4double errMaxNorm, // normalised error
140  G4double hstepCurrent); // current step size
141  // Taking the last step's normalised error, calculate
142  // a step size for the next step.
143  // Limit the next step's size within a range around the current one.
144 
145  inline G4int GetMaxNoSteps() const;
146  inline void SetMaxNoSteps( G4int val);
147  // Modify and Get the Maximum number of Steps that can be
148  // taken for the integration of a single segment -
149  // (ie a single call to AccurateAdvance).
150 
151 //---------------------------------------------------------------------
152 //The following has been introduced by [hackabot] for testing purposes only
153  inline G4int GetTotalNoStepperCalls() const;
154 //---------------------------------------------------------------------
155 
156  public: // without description
157 
158  inline void SetHmin(G4double newval);
159  inline void SetVerboseLevel(G4int newLevel);
160  inline G4double GetVerboseLevel() const;
161 
162  inline G4double GetSmallestFraction() const;
163  void SetSmallestFraction( G4double val );
164 
165 
166 
167 
168  protected: // without description
169  void WarnSmallStepSize( G4double hnext, G4double hstep,
170  G4double h, G4double xDone,
171  G4int noSteps);
172  void WarnTooManyStep( G4double x1start, G4double x2end, G4double xCurrent);
173  void WarnEndPointTooFar (G4double endPointDist,
174  G4double hStepSize ,
175  G4double epsilonRelative,
176  G4int debugFlag);
177  // Issue warnings for undesirable situations
178 
179  void PrintStatus( const G4double* StartArr,
180  G4double xstart,
181  const G4double* CurrentArr,
182  G4double xcurrent,
183  G4double requestStep,
184  G4int subStepNo );
185  void PrintStatus( const G4FieldTrack& StartFT,
186  const G4FieldTrack& CurrentFT,
187  G4double requestStep,
188  G4int subStepNo );
189  void PrintStat_Aux( const G4FieldTrack& aFieldTrack,
190  G4double requestStep,
191  G4double actualStep,
192  G4int subStepNo,
193  G4double subStepSize,
194  G4double dotVelocities );
195  // Verbose output for debugging
196 
197  void PrintStatisticsReport() ;
198  // Report on the number of steps, maximum errors etc.
199 
200 #ifdef QUICK_ADV_TWO
201  G4bool QuickAdvance( G4double yarrin[], // In
202  const G4double dydx[],
203  G4double hstep,
204  G4double yarrout[], // Out
205  G4double& dchord_step, // Out
206  G4double& dyerr ); // in length
207 #endif
208 
209  private:
210 
213  // Private copy constructor and assignment operator.
214 
215  private:
216 
217  // ---------------------------------------------------------------
218  // INVARIANTS
219 
220  G4double fMinimumStep;
221  // Minimum Step allowed in a Step (in absolute units)
222  G4double fSmallestFraction; // Expected range 1e-12 to 5e-15;
223  // Smallest fraction of (existing) curve length - in relative units
224  // below this fraction the current step will be the last
225 
226  const G4int fNoIntegrationVariables; // Number of Variables in integration
227  const G4int fMinNoVars; // Minimum number for FieldTrack
228  const G4int fNoVars; // Full number of variable
229 
230  G4int fMaxNoSteps;
231  static const G4int fMaxStepBase;
232 
233  G4double safety;
234  G4double pshrnk; // exponent for shrinking
235  G4double pgrow; // exponent for growth
236  G4double errcon;
237  // Parameters used to grow and shrink trial stepsize.
238 
239  static const G4double max_stepping_increase;
240  static const G4double max_stepping_decrease;
241  // Maximum stepsize increase/decrease factors.
242 
243  G4int fStatisticsVerboseLevel;
244 
245  // ---------------------------------------------------------------
246  // DEPENDENT Objects
247  G4VFSALIntegrationStepper *pIntStepper;
248 
249  // ---------------------------------------------------------------
250  // STATE
251 
252  G4int fNoTotalSteps, fNoBadSteps, fNoSmallSteps, fNoInitialSmallSteps;
253  G4double fDyerr_max, fDyerr_mx2;
254  G4double fDyerrPos_smTot, fDyerrPos_lgTot, fDyerrVel_lgTot;
255  G4double fSumH_sm, fSumH_lg;
256  // Step Statistics
257 
258  G4int fVerboseLevel; // Verbosity level for printing (debug, ..)
259  // Could be varied during tracking - to help identify issues
260 
261  //For Test Purposes :-
262  G4int TotalNoStepperCalls;
263 
264 };
265 
266 #include "G4FSALIntegrationDriver.icc"
267 
268 #endif /* G4FSALIntegrationDriver_Def */
G4double GetHmin() const
void SetSafety(G4double valS)
void WarnEndPointTooFar(G4double endPointDist, G4double hStepSize, G4double epsilonRelative, G4int debugFlag)
void GetDerivatives(const G4FieldTrack &y_curr, G4double dydx[])
void WarnSmallStepSize(G4double hnext, G4double hstep, G4double h, G4double xDone, G4int noSteps)
void PrintStat_Aux(const G4FieldTrack &aFieldTrack, G4double requestStep, G4double actualStep, G4int subStepNo, G4double subStepSize, G4double dotVelocities)
void SetVerboseLevel(G4int newLevel)
G4double GetSmallestFraction() const
G4double ComputeNewStepSize_WithinLimits(G4double errMaxNorm, G4double hstepCurrent)
G4int GetMaxNoSteps() const
G4double GetVerboseLevel() const
void SetPgrow(G4double valPg)
void SetErrcon(G4double valEc)
static const G4double eps
G4bool AccurateAdvance(G4FieldTrack &y_current, G4double hstep, G4double eps, G4double hinitial=0.0)
int G4int
Definition: G4Types.hh:78
void WarnTooManyStep(G4double x1start, G4double x2end, G4double xCurrent)
G4double GetPshrnk() const
void SetHmin(G4double newval)
G4FSALIntegrationDriver(G4double hminimum, G4VFSALIntegrationStepper *pItsStepper, G4int numberOfComponents=6, G4int statisticsVerbosity=1)
bool G4bool
Definition: G4Types.hh:79
void RenewStepperAndAdjust(G4VFSALIntegrationStepper *pItsStepper)
void SetMaxNoSteps(G4int val)
void SetSmallestFraction(G4double val)
void OneGoodStep(G4double ystart[], G4double dydx[], G4double &x, G4double htry, G4double eps, G4double &hdid, G4double &hnext)
void SetPshrnk(G4double valPs)
G4double ComputeNewStepSize(G4double errMaxNorm, G4double hstepCurrent)
G4int GetTotalNoStepperCalls() const
void ReSetParameters(G4double new_safety=0.9)
G4double Hmin() const
G4double ComputeAndSetErrcon()
G4bool QuickAdvance(G4FieldTrack &y_val, G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr)
double G4double
Definition: G4Types.hh:76
G4double GetPgrow() const
void PrintStatus(const G4double *StartArr, G4double xstart, const G4double *CurrentArr, G4double xcurrent, G4double requestStep, G4int subStepNo)
G4int GetNoTotalSteps() const
G4double GetSafety() const
const G4VFSALIntegrationStepper * GetStepper() const
G4double GetErrcon() const