Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MagInt_Driver Class Reference

G4MagInt_Driver provides a driver that talks to the Integrator Stepper and insures that the error is within acceptable bounds. More...

#include <G4MagIntegratorDriver.hh>

Inheritance diagram for G4MagInt_Driver:

Public Member Functions

 G4MagInt_Driver (G4double hminimum, G4MagIntegratorStepper *pItsStepper, G4int numberOfComponents=6, G4int statisticsVerbosity=0)
 ~G4MagInt_Driver () override
 G4MagInt_Driver (const G4MagInt_Driver &)=delete
G4MagInt_Driveroperator= (const G4MagInt_Driver &)=delete
G4double AdvanceChordLimited (G4FieldTrack &track, G4double stepMax, G4double epsStep, G4double chordDistance) override
void OnStartTracking () override
void OnComputeStep (const G4FieldTrack *=nullptr) override
G4bool DoesReIntegrate () const override
G4bool AccurateAdvance (G4FieldTrack &y_current, G4double hstep, G4double eps, G4double hinitial=0.0) override
G4bool QuickAdvance (G4FieldTrack &y_val, const G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr) override
void StreamInfo (std::ostream &os) const override
G4bool QuickAdvance (G4FieldTrack &y_posvel, const G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr_pos_sq, G4double &dyerr_mom_rel_sq)
G4double GetHmin () const
G4double Hmin () const
G4double GetSafety () const
G4double GetPshrnk () const
G4double GetPgrow () const
G4double GetErrcon () const
void GetDerivatives (const G4FieldTrack &y_curr, G4double dydx[]) const override
void GetDerivatives (const G4FieldTrack &track, G4double dydx[], G4double field[]) const override
G4EquationOfMotionGetEquationOfMotion () override
void SetEquationOfMotion (G4EquationOfMotion *equation) override
void RenewStepperAndAdjust (G4MagIntegratorStepper *pItsStepper) override
void ReSetParameters (G4double new_safety=0.9)
void SetSafety (G4double valS)
void SetPshrnk (G4double valPs)
void SetPgrow (G4double valPg)
void SetErrcon (G4double valEc)
G4double ComputeAndSetErrcon ()
const G4MagIntegratorStepperGetStepper () const override
G4MagIntegratorStepperGetStepper () override
void OneGoodStep (G4double ystart[], const G4double dydx[], G4double &x, G4double htry, G4double eps, G4double &hdid, G4double &hnext)
G4double ComputeNewStepSize (G4double errMaxNorm, G4double hstepCurrent) override
G4double ComputeNewStepSize_WithoutReductionLimit (G4double errMaxNorm, G4double hstepCurrent)
G4double ComputeNewStepSize_WithinLimits (G4double errMaxNorm, G4double hstepCurrent)
G4int GetMaxNoSteps () const
void SetMaxNoSteps (G4int val)
void SetHmin (G4double newval)
void SetVerboseLevel (G4int newLevel) override
G4int GetVerboseLevel () const override
G4double GetSmallestFraction () const
void SetSmallestFraction (G4double val)
Public Member Functions inherited from G4VIntegrationDriver
virtual ~G4VIntegrationDriver ()=default
Public Member Functions inherited from G4ChordFinderDelegate< G4MagInt_Driver >
virtual ~G4ChordFinderDelegate ()
G4double AdvanceChordLimitedImpl (G4FieldTrack &track, G4double hstep, G4double eps, G4double chordDistance)
void ResetStepEstimate ()
G4double GetLastStepEstimateUnc ()
void SetLastStepEstimateUnc (G4double stepEst)
G4int GetNoCalls ()
G4int GetNoTrials ()
G4int GetNoMaxTrials ()
void SetFractions_Last_Next (G4double fractLast=0.90, G4double fractNext=0.95)
void SetFirstFraction (G4double fractFirst)
G4double GetFirstFraction ()
G4double GetFractionLast ()
G4double GetFractionNextEstimate ()
void StreamDelegateInfo (std::ostream &os) const
void TestChordPrint (G4int noTrials, G4int lastStepTrial, G4double dChordStep, G4double fDeltaChord, G4double nextStepTrial)

Protected Member Functions

void WarnSmallStepSize (G4double hnext, G4double hstep, G4double h, G4double xDone, G4int noSteps)
void WarnTooManyStep (G4double x1start, G4double x2end, G4double xCurrent)
void WarnEndPointTooFar (G4double endPointDist, G4double hStepSize, G4double epsilonRelative, G4int debugFlag)
void PrintStatus (const G4double *StartArr, G4double xstart, const G4double *CurrentArr, G4double xcurrent, G4double requestStep, G4int subStepNo)
void PrintStatus (const G4FieldTrack &StartFT, const G4FieldTrack &CurrentFT, G4double requestStep, G4int subStepNo)
void PrintStat_Aux (const G4FieldTrack &aFieldTrack, G4double requestStep, G4double actualStep, G4int subStepNo, G4double subStepSize, G4double dotVelocities)
void PrintStatisticsReport ()

Additional Inherited Members

Static Protected Attributes inherited from G4VIntegrationDriver
static constexpr G4double max_stepping_increase = 5
static constexpr G4double max_stepping_decrease = 0.1

Detailed Description

G4MagInt_Driver provides a driver that talks to the Integrator Stepper and insures that the error is within acceptable bounds.

Definition at line 48 of file G4MagIntegratorDriver.hh.

Constructor & Destructor Documentation

◆ G4MagInt_Driver() [1/2]

G4MagInt_Driver::G4MagInt_Driver ( G4double hminimum,
G4MagIntegratorStepper * pItsStepper,
G4int numberOfComponents = 6,
G4int statisticsVerbosity = 0 )

Constructor for G4MagInt_Driver.

Parameters
[in]hminimumThe minumum allowed step.
[in]pItsStepperPointer to the integrator stepper.
[in]numberOfComponentsThe number of integration variables.
[in]statisticsVerbosityFlag for verbosity.

Definition at line 48 of file G4MagIntegratorDriver.cc.

52 : fNoIntegrationVariables(numComponents),
53 fNoVars( std::max( fNoIntegrationVariables, fMinNoVars )),
54 fStatisticsVerboseLevel(statisticsVerbose)
55{
56 // In order to accomodate "Laboratory Time", which is [7], fMinNoVars=8
57 // is required. For proper time of flight and spin, fMinNoVars must be 12
58
59 RenewStepperAndAdjust( pStepper );
60 fMinimumStep = hminimum;
61
62 fMaxNoSteps = fMaxStepBase / pIntStepper->IntegratorOrder();
63#ifdef G4DEBUG_FIELD
64 fVerboseLevel=2;
65#endif
66
67 if( (fVerboseLevel > 0) || (fStatisticsVerboseLevel > 1) )
68 {
69 G4cout << "MagIntDriver version: Accur-Adv: "
70 << "invE_nS, QuickAdv-2sqrt with Statistics "
71#ifdef G4FLD_STATS
72 << " enabled "
73#else
74 << " disabled "
75#endif
76 << G4endl;
77 }
78}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void RenewStepperAndAdjust(G4MagIntegratorStepper *pItsStepper) override

Referenced by G4MagInt_Driver(), and operator=().

◆ ~G4MagInt_Driver()

G4MagInt_Driver::~G4MagInt_Driver ( )
override

Destructor. Provides statistics if verbosity level is greater than 1.

Definition at line 84 of file G4MagIntegratorDriver.cc.

85{
86 if( fStatisticsVerboseLevel > 1 )
87 {
89 }
90}

◆ G4MagInt_Driver() [2/2]

G4MagInt_Driver::G4MagInt_Driver ( const G4MagInt_Driver & )
delete

Copy constructor and assignment operator not allowed.

Member Function Documentation

◆ AccurateAdvance()

G4bool G4MagInt_Driver::AccurateAdvance ( G4FieldTrack & y_current,
G4double hstep,
G4double eps,
G4double hinitial = 0.0 )
overridevirtual

Advances integration accurately by relative accuracy better than 'eps'.

Parameters
[in,out]y_currentThe current track in field.
[in]hstepProposed step length.
[in]epsRequested accuracy, y_err/hstep.
[in]hinitialInitial minimum integration step.
Returns
true if integration succeeds.

Implements G4VIntegrationDriver.

Definition at line 95 of file G4MagIntegratorDriver.cc.

99{
100 // Runge-Kutta driver with adaptive stepsize control. Integrate starting
101 // values at y_current over hstep x2 with accuracy eps.
102 // On output ystart is replaced by values at the end of the integration
103 // interval. RightHandSide is the right-hand side of ODE system.
104 // The source is similar to odeint routine from NRC p.721-722 .
105
106 G4int nstp, i;
107 G4double x, hnext, hdid, h;
108
109#ifdef G4DEBUG_FIELD
110 G4int no_warnings = 0;
111 static G4int dbg = 1;
112 static G4int nStpPr = 50; // For debug printing of long integrations
113 G4double ySubStepStart[G4FieldTrack::ncompSVEC];
114 G4FieldTrack yFldTrkStart(y_current);
115#endif
116
117 G4double y[G4FieldTrack::ncompSVEC] = {0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.};
118 G4double dydx[G4FieldTrack::ncompSVEC] = {0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.};
119 G4double ystart[G4FieldTrack::ncompSVEC] = {0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.};
120 G4double yEnd[G4FieldTrack::ncompSVEC] = {0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.};
121 G4double x1, x2;
122 G4bool succeeded = true;
123
124 G4double startCurveLength;
125
126 const G4int nvar = fNoVars;
127
128 G4FieldTrack yStartFT(y_current);
129
130 // Ensure that hstep > 0
131 //
132 if( hstep <= 0.0 )
133 {
134 if( hstep == 0.0 )
135 {
136 std::ostringstream message;
137 message << "Proposed step is zero; hstep = " << hstep << " !";
138 G4Exception("G4MagInt_Driver::AccurateAdvance()",
139 "GeomField1001", JustWarning, message);
140 return succeeded;
141 }
142
143 std::ostringstream message;
144 message << "Invalid run condition." << G4endl
145 << "Proposed step is negative; hstep = " << hstep << "." << G4endl
146 << "Requested step cannot be negative! Aborting event.";
147 G4Exception("G4MagInt_Driver::AccurateAdvance()",
148 "GeomField0003", EventMustBeAborted, message);
149 return false;
150 }
151
152 y_current.DumpToArray( ystart );
153
154 startCurveLength= y_current.GetCurveLength();
155 x1= startCurveLength;
156 x2= x1 + hstep;
157
158 if ( (hinitial > 0.0) && (hinitial < hstep)
159 && (hinitial > perMillion * hstep) )
160 {
161 h = hinitial;
162 }
163 else // Initial Step size "h" defaults to the full interval
164 {
165 h = hstep;
166 }
167
168 x = x1;
169
170 for ( i=0; i<nvar; ++i) { y[i] = ystart[i]; }
171
172 G4bool lastStep= false;
173 nstp = 1;
174
175 do
176 {
177 G4ThreeVector StartPos( y[0], y[1], y[2] );
178
179#ifdef G4DEBUG_FIELD
180 G4double xSubStepStart= x;
181 for (i=0; i<nvar; ++i) { ySubStepStart[i] = y[i]; }
182 yFldTrkStart.LoadFromArray(y, fNoIntegrationVariables);
183 yFldTrkStart.SetCurveLength(x);
184#endif
185
186 pIntStepper->RightHandSide( y, dydx );
187 ++fNoTotalSteps;
188
189 // Perform the Integration
190 //
191 if( h > fMinimumStep )
192 {
193 OneGoodStep(y,dydx,x,h,eps,hdid,hnext) ;
194 //--------------------------------------
195#ifdef G4DEBUG_FIELD
196 if (dbg>2)
197 {
198 // PrintStatus( ySubStepStart, xSubStepStart, y, x, h, nstp); // Only
199 G4DriverReporter::PrintStatus( ySubStepStart, xSubStepStart, y, x, h, nstp, nvar);
200 }
201#endif
202 }
203 else
204 {
205 G4FieldTrack yFldTrk( G4ThreeVector(0,0,0),
206 G4ThreeVector(0,0,0), 0., 0., 0., 0. );
207 G4double dchord_step, dyerr, dyerr_len; // What to do with these ?
208 yFldTrk.LoadFromArray(y, fNoIntegrationVariables);
209 yFldTrk.SetCurveLength( x );
210
211 QuickAdvance( yFldTrk, dydx, h, dchord_step, dyerr_len );
212 //-----------------------------------------------------
213
214 yFldTrk.DumpToArray(y);
215
216#ifdef G4FLD_STATS
217 ++fNoSmallSteps;
218 if ( dyerr_len > fDyerr_max ) { fDyerr_max = dyerr_len; }
219 fDyerrPos_smTot += dyerr_len;
220 fSumH_sm += h; // Length total for 'small' steps
221 if (nstp==1) { ++fNoInitialSmallSteps; }
222#endif
223#ifdef G4DEBUG_FIELD
224 if (dbg>1)
225 {
226 if(fNoSmallSteps<2) { PrintStatus(ySubStepStart, x1, y, x, h, -nstp); }
227 G4cout << "Another sub-min step, no " << fNoSmallSteps
228 << " of " << fNoTotalSteps << " this time " << nstp << G4endl;
229 PrintStatus( ySubStepStart, x1, y, x, h, nstp); // Only this
230 G4cout << " dyerr= " << dyerr_len << " relative = " << dyerr_len / h
231 << " epsilon= " << eps << " hstep= " << hstep
232 << " h= " << h << " hmin= " << fMinimumStep << G4endl;
233 }
234#endif
235 if( h == 0.0 )
236 {
237 G4Exception("G4MagInt_Driver::AccurateAdvance()",
238 "GeomField0003", FatalException,
239 "Integration Step became Zero!");
240 }
241 dyerr = dyerr_len / h;
242 hdid = h;
243 x += hdid;
244
245 // Compute suggested new step
246 hnext = ComputeNewStepSize( dyerr/eps, h);
247 }
248
249 G4ThreeVector EndPos( y[0], y[1], y[2] );
250
251#ifdef G4DEBUG_FIELD
252 if( (dbg>0) && (dbg<=2) && (nstp>nStpPr))
253 {
254 if( nstp==nStpPr ) { G4cout << "***** Many steps ****" << G4endl; }
255 G4cout << "MagIntDrv: " ;
256 G4cout << "hdid=" << std::setw(12) << hdid << " "
257 << "hnext=" << std::setw(12) << hnext << " "
258 << "hstep=" << std::setw(12) << hstep << " (requested) "
259 << G4endl;
260 PrintStatus( ystart, x1, y, x, h, (nstp==nStpPr) ? -nstp: nstp);
261 }
262#endif
263
264 // Check the endpoint
265 G4double endPointDist= (EndPos-StartPos).mag();
266 if ( endPointDist >= hdid*(1.+perMillion) )
267 {
268 ++fNoBadSteps;
269
270 // Issue a warning only for gross differences -
271 // we understand how small difference occur.
272 if ( endPointDist >= hdid*(1.+perThousand) )
273 {
274#ifdef G4DEBUG_FIELD
275 if (dbg)
276 {
277 WarnEndPointTooFar ( endPointDist, hdid, eps, dbg );
278 G4cerr << " Total steps: bad " << fNoBadSteps
279 << " current h= " << hdid << G4endl;
280 PrintStatus( ystart, x1, y, x, hstep, no_warnings?nstp:-nstp);
281 }
282 ++no_warnings;
283#endif
284 }
285 }
286
287 // Avoid numerous small last steps
288 if( (h < eps * hstep) || (h < fSmallestFraction * startCurveLength) )
289 {
290 // No more integration -- the next step will not happen
291 lastStep = true;
292 }
293 else
294 {
295 // Check the proposed next stepsize
296 if(std::fabs(hnext) <= Hmin())
297 {
298#ifdef G4DEBUG_FIELD
299 // If simply a very small interval is being integrated, do not warn
300 if( (x < x2 * (1-eps) ) && // The last step can be small: OK
301 (std::fabs(hstep) > Hmin()) ) // and if we are asked, it's OK
302 {
303 if(dbg>0)
304 {
305 WarnSmallStepSize( hnext, hstep, h, x-x1, nstp );
306 PrintStatus( ystart, x1, y, x, hstep, no_warnings?nstp:-nstp);
307 }
308 ++no_warnings;
309 }
310#endif
311 // Make sure that the next step is at least Hmin.
312 h = Hmin();
313 }
314 else
315 {
316 h = hnext;
317 }
318
319 // Ensure that the next step does not overshoot
320 if ( x+h > x2 )
321 { // When stepsize overshoots, decrease it!
322 h = x2 - x ; // Must cope with difficult rounding-error
323 } // issues if hstep << x2
324
325 if ( h == 0.0 )
326 {
327 // Cannot progress - accept this as last step - by default
328 lastStep = true;
329#ifdef G4DEBUG_FIELD
330 if (dbg>2)
331 {
332 int prec= G4cout.precision(12);
333 G4cout << "Warning: G4MagIntegratorDriver::AccurateAdvance"
334 << G4endl
335 << " Integration step 'h' became "
336 << h << " due to roundoff. " << G4endl
337 << " Calculated as difference of x2= "<< x2 << " and x=" << x
338 << " Forcing termination of advance." << G4endl;
339 G4cout.precision(prec);
340 }
341#endif
342 }
343 }
344 } while ( ((++nstp)<=fMaxNoSteps) && (x < x2) && (!lastStep) );
345 // Loop checking, 07.10.2016, J. Apostolakis
346
347 // Have we reached the end ?
348 // --> a better test might be x-x2 > an_epsilon
349
350 succeeded = (x>=x2); // If it was a "forced" last step
351
352 for (i=0; i<nvar; ++i) { yEnd[i] = y[i]; }
353
354 // Put back the values.
355 y_current.LoadFromArray( yEnd, fNoIntegrationVariables );
356 y_current.SetCurveLength( x );
357
358 if(nstp > fMaxNoSteps)
359 {
360 succeeded = false;
361#ifdef G4DEBUG_FIELD
362 ++no_warnings;
363 if (dbg)
364 {
365 WarnTooManyStep( x1, x2, x ); // Issue WARNING
366 PrintStatus( yEnd, x1, y, x, hstep, -nstp);
367 }
368#endif
369 }
370
371#ifdef G4DEBUG_FIELD
372 if( dbg && no_warnings )
373 {
374 G4cerr << "G4MagIntegratorDriver exit status: no-steps " << nstp << G4endl;
375 PrintStatus( yEnd, x1, y, x, hstep, nstp);
376 }
377#endif
378
379 return succeeded;
380} // end of AccurateAdvance ...........................
@ JustWarning
@ FatalException
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
static void PrintStatus(const G4double *StartArr, G4double xstart, const G4double *CurrentArr, G4double xcurrent, G4double requestStep, unsigned int subStepNo, unsigned int noIntegrationVariables)
G4double GetCurveLength() const
void SetCurveLength(G4double nCurve_s)
void DumpToArray(G4double valArr[ncompSVEC]) const
void LoadFromArray(const G4double valArr[ncompSVEC], G4int noVarsIntegrated)
void PrintStatus(const G4double *StartArr, G4double xstart, const G4double *CurrentArr, G4double xcurrent, G4double requestStep, G4int subStepNo)
G4bool QuickAdvance(G4FieldTrack &y_val, const G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr) override
G4double Hmin() const
void WarnSmallStepSize(G4double hnext, G4double hstep, G4double h, G4double xDone, G4int noSteps)
void WarnEndPointTooFar(G4double endPointDist, G4double hStepSize, G4double epsilonRelative, G4int debugFlag)
void OneGoodStep(G4double ystart[], const G4double dydx[], G4double &x, G4double htry, G4double eps, G4double &hdid, G4double &hnext)
G4double ComputeNewStepSize(G4double errMaxNorm, G4double hstepCurrent) override
void WarnTooManyStep(G4double x1start, G4double x2end, G4double xCurrent)

◆ AdvanceChordLimited()

G4double G4MagInt_Driver::AdvanceChordLimited ( G4FieldTrack & track,
G4double stepMax,
G4double epsStep,
G4double chordDistance )
inlineoverridevirtual

Computes the step to take, based on chord limits.

Parameters
[in,out]trackThe current track in field.
[in]stepMaxProposed maximum step length.
[in]epsStepRequested accuracy, y_err/hstep.
[in]chordDistanceMaximum sagitta distance.
Returns
The length of step taken.

Implements G4VIntegrationDriver.

◆ ComputeAndSetErrcon()

G4double G4MagInt_Driver::ComputeAndSetErrcon ( )
inline

◆ ComputeNewStepSize()

G4double G4MagInt_Driver::ComputeNewStepSize ( G4double errMaxNorm,
G4double hstepCurrent )
overridevirtual

Takes the last step's normalised error and calculates a step size for the next step. Does it limit the next step's size within a factor of the current? – DOES NOT limit for very bad steps – DOES limit for very good (x5).

Implements G4VIntegrationDriver.

Definition at line 735 of file G4MagIntegratorDriver.cc.

738{
739 // Legacy behaviour:
740 return ComputeNewStepSize_WithoutReductionLimit( errMaxNorm, hstepCurrent );
741 // 'Improved' behaviour - at least more consistent with other step estimates:
742 // return ComputeNewStepSize_WithinLimits( errMaxNorm, hstepCurrent );
743}
G4double ComputeNewStepSize_WithoutReductionLimit(G4double errMaxNorm, G4double hstepCurrent)

Referenced by AccurateAdvance(), and QuickAdvance().

◆ ComputeNewStepSize_WithinLimits()

G4double G4MagInt_Driver::ComputeNewStepSize_WithinLimits ( G4double errMaxNorm,
G4double hstepCurrent )

Taking the last step's normalised error, calculates a step size for the next step. Limits the next step's size within a range around the current one.

Definition at line 751 of file G4MagIntegratorDriver.cc.

754{
755 G4double hnew;
756
757 // Compute size of next Step for a failed step
758 if (errMaxNorm > 1.0 )
759 {
760 // Step failed; compute the size of retrial Step.
761 hnew = GetSafety()*hstepCurrent*std::pow(errMaxNorm,GetPshrnk()) ;
762
763 if (hnew < max_stepping_decrease*hstepCurrent)
764 {
765 hnew = max_stepping_decrease*hstepCurrent ;
766 // reduce stepsize, but no more
767 // than this factor (value= 1/10)
768 }
769 }
770 else
771 {
772 // Compute size of next Step for a successful step
773 if (errMaxNorm > errcon)
774 { hnew = GetSafety()*hstepCurrent*std::pow(errMaxNorm,GetPgrow()); }
775 else // No more than a factor of 5 increase
776 { hnew = max_stepping_increase * hstepCurrent; }
777 }
778 return hnew;
779}
G4double GetPshrnk() const
G4double GetSafety() const
G4double GetPgrow() const
static constexpr G4double max_stepping_decrease
static constexpr G4double max_stepping_increase

◆ ComputeNewStepSize_WithoutReductionLimit()

G4double G4MagInt_Driver::ComputeNewStepSize_WithoutReductionLimit ( G4double errMaxNorm,
G4double hstepCurrent )

Taking the last step's normalised error, calculates a step size for the next step. Does not limit the next step's size within a factor of the current one when reducing the size, i.e. for badly failing steps.

Definition at line 706 of file G4MagIntegratorDriver.cc.

709{
710 G4double hnew;
711
712 // Compute size of next Step for a failed step
713 if(errMaxNorm > 1.0 )
714 {
715 // Step failed; compute the size of retrial Step.
716 hnew = GetSafety()*hstepCurrent*std::pow(errMaxNorm,GetPshrnk()) ;
717 }
718 else if(errMaxNorm > 0.0 )
719 {
720 // Compute size of next Step for a successful step
721 hnew = GetSafety()*hstepCurrent*std::pow(errMaxNorm,GetPgrow()) ;
722 }
723 else
724 {
725 // if error estimate is zero (possible) or negative (dubious)
726 hnew = max_stepping_increase * hstepCurrent;
727 }
728
729 return hnew;
730}

Referenced by ComputeNewStepSize().

◆ DoesReIntegrate()

G4bool G4MagInt_Driver::DoesReIntegrate ( ) const
inlineoverridevirtual

The driver implements re-integration, so returns true.

Implements G4VIntegrationDriver.

Definition at line 102 of file G4MagIntegratorDriver.hh.

102{ return true; }

Referenced by PrintInfo(), and StreamInfo().

◆ GetDerivatives() [1/2]

void G4MagInt_Driver::GetDerivatives ( const G4FieldTrack & track,
G4double dydx[],
G4double field[] ) const
overridevirtual

Implements G4VIntegrationDriver.

Definition at line 980 of file G4MagIntegratorDriver.cc.

983{
985 track.DumpToArray(ytemp);
986 pIntStepper->RightHandSide(ytemp, dydx, field);
987}

◆ GetDerivatives() [2/2]

void G4MagInt_Driver::GetDerivatives ( const G4FieldTrack & y_curr,
G4double dydx[] ) const
overridevirtual

Implements G4VIntegrationDriver.

◆ GetEquationOfMotion()

G4EquationOfMotion * G4MagInt_Driver::GetEquationOfMotion ( )
overridevirtual

Getter and setter for the equation of motion.

Implements G4VIntegrationDriver.

Definition at line 989 of file G4MagIntegratorDriver.cc.

990{
991 return pIntStepper->GetEquationOfMotion();
992}

◆ GetErrcon()

G4double G4MagInt_Driver::GetErrcon ( ) const
inline

Referenced by PrintInfo().

◆ GetHmin()

G4double G4MagInt_Driver::GetHmin ( ) const
inline

Accessors.

Referenced by PrintInfo().

◆ GetMaxNoSteps()

G4int G4MagInt_Driver::GetMaxNoSteps ( ) const
inline

Modifier and accessor for the maximum number of steps that can be taken for the integration of a single segment, i.e. a single call to AccurateAdvance().

Referenced by PrintInfo().

◆ GetPgrow()

◆ GetPshrnk()

G4double G4MagInt_Driver::GetPshrnk ( ) const
inline

◆ GetSafety()

G4double G4MagInt_Driver::GetSafety ( ) const
inline

◆ GetSmallestFraction()

G4double G4MagInt_Driver::GetSmallestFraction ( ) const
inline

Referenced by PrintInfo().

◆ GetStepper() [1/2]

const G4MagIntegratorStepper * G4MagInt_Driver::GetStepper ( ) const
overridevirtual

Accessors for the integrator stepper.

Implements G4VIntegrationDriver.

Definition at line 999 of file G4MagIntegratorDriver.cc.

1000{
1001 return pIntStepper;
1002}

◆ GetStepper() [2/2]

G4MagIntegratorStepper * G4MagInt_Driver::GetStepper ( )
overridevirtual

Implements G4VIntegrationDriver.

Definition at line 1004 of file G4MagIntegratorDriver.cc.

1005{
1006 return pIntStepper;
1007}

◆ GetVerboseLevel()

G4int G4MagInt_Driver::GetVerboseLevel ( ) const
overridevirtual

Implements G4VIntegrationDriver.

Referenced by PrintInfo().

◆ Hmin()

G4double G4MagInt_Driver::Hmin ( ) const
inline

◆ OnComputeStep()

void G4MagInt_Driver::OnComputeStep ( const G4FieldTrack * = nullptr)
inlineoverridevirtual

Dispatch interface method for computing step. Does nothing here.

Implements G4VIntegrationDriver.

Definition at line 97 of file G4MagIntegratorDriver.hh.

97{}

◆ OneGoodStep()

void G4MagInt_Driver::OneGoodStep ( G4double ystart[],
const G4double dydx[],
G4double & x,
G4double htry,
G4double eps,
G4double & hdid,
G4double & hnext )

Takes one Step that is as large as possible while satisfying the accuracy criterion of: yerr < eps * |y_end-y_start|.

Parameters
[in,out]ystartThe current track state, y.
[in]dydxThe derivatives array.
[in,out]xStep start, x.
[in]htryStep to attempt.
[in]epsThe relative accuracy.
[out]hdidStep achieved.
[out]hnextProposed next step.
Returns
true if integration succeeds.

Definition at line 470 of file G4MagIntegratorDriver.cc.

491{
492 G4double errmax_sq;
493 G4double h, htemp, xnew ;
494
496
497 h = htry ; // Set stepsize to the initial trial value
498
499 G4double inv_eps_vel_sq = 1.0 / (eps_rel_max*eps_rel_max);
500
501 G4double errpos_sq = 0.0; // square of displacement error
502 G4double errvel_sq = 0.0; // square of momentum vector difference
503 G4double errspin_sq = 0.0; // square of spin vector difference
504
505 const G4int max_trials=100;
506
507 G4ThreeVector Spin(y[9],y[10],y[11]);
508 G4double spin_mag2 = Spin.mag2();
509 G4bool hasSpin = (spin_mag2 > 0.0);
510
511 for (G4int iter=0; iter<max_trials; ++iter)
512 {
513 pIntStepper-> Stepper(y,dydx,h,ytemp,yerr);
514 // *******
515 G4double eps_pos = eps_rel_max * std::max(h, fMinimumStep);
516 G4double inv_eps_pos_sq = 1.0 / (eps_pos*eps_pos);
517
518 // Evaluate accuracy
519 //
520 errpos_sq = sqr(yerr[0]) + sqr(yerr[1]) + sqr(yerr[2]) ;
521 errpos_sq *= inv_eps_pos_sq; // Scale relative to required tolerance
522
523 // Accuracy for momentum
524 G4double magvel_sq= sqr(y[3]) + sqr(y[4]) + sqr(y[5]) ;
525 G4double sumerr_sq = sqr(yerr[3]) + sqr(yerr[4]) + sqr(yerr[5]) ;
526 if( magvel_sq > 0.0 )
527 {
528 errvel_sq = sumerr_sq / magvel_sq;
529 }
530 else
531 {
532 std::ostringstream message;
533 message << "Found case of zero momentum." << G4endl
534 << "- iteration= " << iter << "; h= " << h;
535 G4Exception("G4MagInt_Driver::OneGoodStep()",
536 "GeomField1001", JustWarning, message);
537 errvel_sq = sumerr_sq;
538 }
539 errvel_sq *= inv_eps_vel_sq;
540 errmax_sq = std::max( errpos_sq, errvel_sq ); // Square of maximum error
541
542 if( hasSpin )
543 {
544 // Accuracy for spin
545 errspin_sq = ( sqr(yerr[9]) + sqr(yerr[10]) + sqr(yerr[11]) )
546 / spin_mag2; // ( sqr(y[9]) + sqr(y[10]) + sqr(y[11]) );
547 errspin_sq *= inv_eps_vel_sq;
548 errmax_sq = std::max( errmax_sq, errspin_sq );
549 }
550
551 if ( errmax_sq <= 1.0 ) { break; } // Step succeeded.
552
553 // Step failed; compute the size of retrial Step.
554 htemp = GetSafety() * h * std::pow( errmax_sq, 0.5*GetPshrnk() );
555
556 if (htemp >= 0.1*h) { h = htemp; } // Truncation error too large,
557 else { h = 0.1*h; } // reduce stepsize, but no more
558 // than a factor of 10
559 xnew = x + h;
560 if(xnew == x)
561 {
562 std::ostringstream message;
563 message << "Stepsize underflow in Stepper !" << G4endl
564 << "- Step's start x=" << x << " and end x= " << xnew
565 << " are equal !! " << G4endl
566 << " Due to step-size= " << h
567 << ". Note that input step was " << htry;
568 G4Exception("G4MagInt_Driver::OneGoodStep()",
569 "GeomField1001", JustWarning, message);
570 break;
571 }
572 }
573
574 // Compute size of next Step
575 if (errmax_sq > errcon*errcon)
576 {
577 hnext = GetSafety()*h*std::pow(errmax_sq, 0.5*GetPgrow());
578 }
579 else
580 {
581 hnext = max_stepping_increase*h ; // No more than a factor of 5 increase
582 }
583 x += (hdid = h);
584
585 for(G4int k=0; k<fNoIntegrationVariables; ++k) { y[k] = ytemp[k]; }
586
587 return;
588}
T sqr(const T &x)
Definition templates.hh:128

Referenced by AccurateAdvance().

◆ OnStartTracking()

void G4MagInt_Driver::OnStartTracking ( )
inlineoverridevirtual

Dispatch interface method for initialisation/reset of driver.

Implements G4VIntegrationDriver.

◆ operator=()

G4MagInt_Driver & G4MagInt_Driver::operator= ( const G4MagInt_Driver & )
delete

◆ PrintStat_Aux()

void G4MagInt_Driver::PrintStat_Aux ( const G4FieldTrack & aFieldTrack,
G4double requestStep,
G4double actualStep,
G4int subStepNo,
G4double subStepSize,
G4double dotVelocities )
protected

Definition at line 871 of file G4MagIntegratorDriver.cc.

877{
878 const G4ThreeVector Position = aFieldTrack.GetPosition();
879 const G4ThreeVector UnitVelocity = aFieldTrack.GetMomentumDir();
880
881 if( subStepNo >= 0)
882 {
883 G4cout << std::setw( 5) << subStepNo << " ";
884 }
885 else
886 {
887 G4cout << std::setw( 5) << "Start" << " ";
888 }
889 G4double curveLen= aFieldTrack.GetCurveLength();
890 G4cout << std::setw( 7) << curveLen;
891 G4cout << std::setw( 9) << Position.x() << " "
892 << std::setw( 9) << Position.y() << " "
893 << std::setw( 9) << Position.z() << " "
894 << std::setw( 8) << UnitVelocity.x() << " "
895 << std::setw( 8) << UnitVelocity.y() << " "
896 << std::setw( 8) << UnitVelocity.z() << " ";
897 G4long oldprec= G4cout.precision(3);
898 G4cout << std::setw( 8) << UnitVelocity.mag2()-1.0 << " ";
899 G4cout.precision(6);
900 G4cout << std::setw(10) << dotVeloc_StartCurr << " ";
901 G4cout.precision(oldprec);
902 G4cout << std::setw( 7) << aFieldTrack.GetKineticEnergy();
903 G4cout << std::setw(12) << step_len << " ";
904
905 static G4ThreadLocal G4double oldCurveLength = 0.0;
906 static G4ThreadLocal G4double oldSubStepLength = 0.0;
907 static G4ThreadLocal G4int oldSubStepNo = -1;
908
909 G4double subStep_len = 0.0;
910 if( curveLen > oldCurveLength )
911 {
912 subStep_len= curveLen - oldCurveLength;
913 }
914 else if (subStepNo == oldSubStepNo)
915 {
916 subStep_len= oldSubStepLength;
917 }
918 oldCurveLength= curveLen;
919 oldSubStepLength= subStep_len;
920
921 G4cout << std::setw(12) << subStep_len << " ";
922 G4cout << std::setw(12) << subStepSize << " ";
923 if( requestStep != -1.0 )
924 {
925 G4cout << std::setw( 9) << requestStep << " ";
926 }
927 else
928 {
929 G4cout << std::setw( 9) << " InitialStep " << " ";
930 }
931 G4cout << G4endl;
932}
long G4long
Definition G4Types.hh:87
double z() const
double x() const
double mag2() const
double y() const
const G4ThreeVector & GetMomentumDir() const
G4double GetKineticEnergy() const
G4ThreeVector GetPosition() const
#define G4ThreadLocal
Definition tls.hh:77

Referenced by PrintStatus().

◆ PrintStatisticsReport()

void G4MagInt_Driver::PrintStatisticsReport ( )
protected

Reports on the number of steps, maximum errors etc.

Definition at line 936 of file G4MagIntegratorDriver.cc.

937{
938 G4int noPrecBig = 6;
939 G4long oldPrec = G4cout.precision(noPrecBig);
940
941 G4cout << "G4MagInt_Driver Statistics of steps undertaken. " << G4endl;
942 G4cout << "G4MagInt_Driver: Number of Steps: "
943 << " Total= " << fNoTotalSteps
944 << " Bad= " << fNoBadSteps
945 << " Small= " << fNoSmallSteps
946 << " Non-initial small= " << (fNoSmallSteps-fNoInitialSmallSteps)
947 << G4endl;
948 G4cout.precision(oldPrec);
949}

Referenced by ~G4MagInt_Driver().

◆ PrintStatus() [1/2]

void G4MagInt_Driver::PrintStatus ( const G4double * StartArr,
G4double xstart,
const G4double * CurrentArr,
G4double xcurrent,
G4double requestStep,
G4int subStepNo )
protected

Loggers for verbosity printouts.

Definition at line 783 of file G4MagIntegratorDriver.cc.

793{
794 G4FieldTrack StartFT(G4ThreeVector(0,0,0),
795 G4ThreeVector(0,0,0), 0., 0., 0., 0. );
796 G4FieldTrack CurrentFT (StartFT);
797
798 StartFT.LoadFromArray( StartArr, fNoIntegrationVariables);
799 StartFT.SetCurveLength( xstart);
800 CurrentFT.LoadFromArray( CurrentArr, fNoIntegrationVariables);
801 CurrentFT.SetCurveLength( xcurrent );
802
803 PrintStatus(StartFT, CurrentFT, requestStep, subStepNo );
804}

Referenced by AccurateAdvance(), PrintStatus(), and QuickAdvance().

◆ PrintStatus() [2/2]

void G4MagInt_Driver::PrintStatus ( const G4FieldTrack & StartFT,
const G4FieldTrack & CurrentFT,
G4double requestStep,
G4int subStepNo )
protected

Definition at line 808 of file G4MagIntegratorDriver.cc.

812{
813 G4int verboseLevel= fVerboseLevel;
814 const G4int noPrecision = 5;
815 G4long oldPrec= G4cout.precision(noPrecision);
816 // G4cout.setf(ios_base::fixed,ios_base::floatfield);
817
818 const G4ThreeVector StartPosition= StartFT.GetPosition();
819 const G4ThreeVector StartUnitVelocity= StartFT.GetMomentumDir();
820 const G4ThreeVector CurrentPosition= CurrentFT.GetPosition();
821 const G4ThreeVector CurrentUnitVelocity= CurrentFT.GetMomentumDir();
822
823 G4double DotStartCurrentVeloc= StartUnitVelocity.dot(CurrentUnitVelocity);
824
825 G4double step_len= CurrentFT.GetCurveLength() - StartFT.GetCurveLength();
826 G4double subStepSize = step_len;
827
828 if( (subStepNo <= 1) || (verboseLevel > 3) )
829 {
830 subStepNo = - subStepNo; // To allow printing banner
831
832 G4cout << std::setw( 6) << " " << std::setw( 25)
833 << " G4MagInt_Driver: Current Position and Direction" << " "
834 << G4endl;
835 G4cout << std::setw( 5) << "Step#" << " "
836 << std::setw( 7) << "s-curve" << " "
837 << std::setw( 9) << "X(mm)" << " "
838 << std::setw( 9) << "Y(mm)" << " "
839 << std::setw( 9) << "Z(mm)" << " "
840 << std::setw( 8) << " N_x " << " "
841 << std::setw( 8) << " N_y " << " "
842 << std::setw( 8) << " N_z " << " "
843 << std::setw( 8) << " N^2-1 " << " "
844 << std::setw(10) << " N(0).N " << " "
845 << std::setw( 7) << "KinEner " << " "
846 << std::setw(12) << "Track-l" << " " // Add the Sub-step ??
847 << std::setw(12) << "Step-len" << " "
848 << std::setw(12) << "Step-len" << " "
849 << std::setw( 9) << "ReqStep" << " "
850 << G4endl;
851 }
852
853 if( (subStepNo <= 0) )
854 {
855 PrintStat_Aux( StartFT, requestStep, 0.,
856 0, 0.0, 1.0);
857 }
858
859 if( verboseLevel <= 3 )
860 {
861 G4cout.precision(noPrecision);
862 PrintStat_Aux( CurrentFT, requestStep, step_len,
863 subStepNo, subStepSize, DotStartCurrentVeloc );
864 }
865
866 G4cout.precision(oldPrec);
867}
double dot(const Hep3Vector &) const
void PrintStat_Aux(const G4FieldTrack &aFieldTrack, G4double requestStep, G4double actualStep, G4int subStepNo, G4double subStepSize, G4double dotVelocities)

◆ QuickAdvance() [1/2]

G4bool G4MagInt_Driver::QuickAdvance ( G4FieldTrack & y_posvel,
const G4double dydx[],
G4double hstep,
G4double & dchord_step,
G4double & dyerr_pos_sq,
G4double & dyerr_mom_rel_sq )

Attempts one integration step, and returns estimated error 'dyerr'. It does not ensure accuracy.

Parameters
[in,out]y_posvelThe current track in field.
[in]dydxdydx array.
[in]hstepProposed step length.
[out]dchord_stepEstimated sagitta distance.
[out]dyerr_pos_sqEstimated error in position.
[out]dyerr_mom_rel_sqEstimated error in momentum (normalised: Delta_Integration(p^2)/(p^2)).
Returns
true if integration succeeds.

Definition at line 594 of file G4MagIntegratorDriver.cc.

600{
601 G4Exception("G4MagInt_Driver::QuickAdvance()", "GeomField0001",
602 FatalException, "Not yet implemented.");
603
604 // Use the parameters of this method, to please compiler
605 //
606 dchord_step = dyerr_pos_sq = hstep * hstep * dydx[0];
607 dyerr_mom_rel_sq = y_posvel.GetPosition().mag2();
608 return true;
609}

◆ QuickAdvance() [2/2]

G4bool G4MagInt_Driver::QuickAdvance ( G4FieldTrack & y_val,
const G4double dydx[],
G4double hstep,
G4double & dchord_step,
G4double & dyerr )
overridevirtual

Attempts one integration step, and returns estimated error 'dyerr'. It does not ensure accuracy.

Parameters
[in,out]y_valThe current track in field.
[in]dydxdydx array.
[in]hstepProposed step length.
[out]dchord_stepEstimated sagitta distance.
[out]dyerrEstimated error.
Returns
true if integration succeeds.

Reimplemented from G4VIntegrationDriver.

Definition at line 613 of file G4MagIntegratorDriver.cc.

618{
619 G4double dyerr_pos_sq, dyerr_mom_rel_sq;
622 G4double s_start;
623 G4double dyerr_mom_sq, vel_mag_sq, inv_vel_mag_sq;
624
625 // Move data into array
626 y_posvel.DumpToArray( yarrin ); // yarrin <== y_posvel
627 s_start = y_posvel.GetCurveLength();
628
629 // Do an Integration Step
630 pIntStepper-> Stepper(yarrin, dydx, hstep, yarrout, yerr_vec) ;
631
632 // Estimate curve-chord distance
633 dchord_step= pIntStepper-> DistChord();
634
635 // Put back the values. yarrout ==> y_posvel
636 y_posvel.LoadFromArray( yarrout, fNoIntegrationVariables );
637 y_posvel.SetCurveLength( s_start + hstep );
638
639#ifdef G4DEBUG_FIELD
640 if(fVerboseLevel>2)
641 {
642 G4cout << "G4MagIntDrv: Quick Advance" << G4endl;
643 PrintStatus( yarrin, s_start, yarrout, s_start+hstep, hstep, 1);
644 }
645#endif
646
647 // A single measure of the error
648 // TO-DO : account for energy, spin, ... ?
649 vel_mag_sq = ( sqr(yarrout[3])+sqr(yarrout[4])+sqr(yarrout[5]) );
650 inv_vel_mag_sq = 1.0 / vel_mag_sq;
651 dyerr_pos_sq = ( sqr(yerr_vec[0])+sqr(yerr_vec[1])+sqr(yerr_vec[2]));
652 dyerr_mom_sq = ( sqr(yerr_vec[3])+sqr(yerr_vec[4])+sqr(yerr_vec[5]));
653 dyerr_mom_rel_sq = dyerr_mom_sq * inv_vel_mag_sq;
654
655 // Calculate also the change in the momentum squared also ???
656 // G4double veloc_square = y_posvel.GetVelocity().mag2();
657 // ...
658
659#ifdef RETURN_A_NEW_STEP_LENGTH
660 // The following step cannot be done here because "eps" is not known.
661 dyerr_len = std::sqrt( dyerr_len_sq );
662 dyerr_len_sq /= eps ;
663
664 // Look at the velocity deviation ?
665 // sqr(yerr_vec[3])+sqr(yerr_vec[4])+sqr(yerr_vec[5]));
666
667 // Set suggested new step
668 hstep = ComputeNewStepSize( dyerr_len, hstep);
669#endif
670
671 if( dyerr_pos_sq > ( dyerr_mom_rel_sq * sqr(hstep) ) )
672 {
673 dyerr = std::sqrt(dyerr_pos_sq);
674 }
675 else
676 {
677 // Scale it to the current step size - for now
678 dyerr = std::sqrt(dyerr_mom_rel_sq) * hstep;
679 }
680
681 return true;
682}

Referenced by AccurateAdvance().

◆ RenewStepperAndAdjust()

void G4MagInt_Driver::RenewStepperAndAdjust ( G4MagIntegratorStepper * pItsStepper)
overridevirtual

Sets a new stepper 'pItsStepper' for this driver. Then it calls ResetParameters() to update its parameters accordingly.

Reimplemented from G4VIntegrationDriver.

Definition at line 1009 of file G4MagIntegratorDriver.cc.

1011{
1012 pIntStepper = pItsStepper;
1014}
void ReSetParameters(G4double new_safety=0.9)

Referenced by G4MagInt_Driver().

◆ ReSetParameters()

void G4MagInt_Driver::ReSetParameters ( G4double new_safety = 0.9)
inline

Resets the qarameters according to the new provided safety value. i) sets the exponents (pgrow & pshrnk), using the current order; ii) sets the safety and calculates "errcon" according to the above values.

Referenced by RenewStepperAndAdjust().

◆ SetEquationOfMotion()

void G4MagInt_Driver::SetEquationOfMotion ( G4EquationOfMotion * equation)
overridevirtual

Setter and getter for the equation of motion.

Implements G4VIntegrationDriver.

Definition at line 994 of file G4MagIntegratorDriver.cc.

995{
996 pIntStepper->SetEquationOfMotion(equation);
997}

◆ SetErrcon()

void G4MagInt_Driver::SetErrcon ( G4double valEc)
inline

◆ SetHmin()

void G4MagInt_Driver::SetHmin ( G4double newval)
inline

More modifiers and accessors.

◆ SetMaxNoSteps()

void G4MagInt_Driver::SetMaxNoSteps ( G4int val)
inline

◆ SetPgrow()

void G4MagInt_Driver::SetPgrow ( G4double valPg)
inline

◆ SetPshrnk()

void G4MagInt_Driver::SetPshrnk ( G4double valPs)
inline

◆ SetSafety()

void G4MagInt_Driver::SetSafety ( G4double valS)
inline

Modifiers. When setting safety or pgrow, errcon will be set to a compatible value.

◆ SetSmallestFraction()

void G4MagInt_Driver::SetSmallestFraction ( G4double val)

Definition at line 953 of file G4MagIntegratorDriver.cc.

954{
955 if( (newFraction > 1.e-16) && (newFraction < 1e-8) )
956 {
957 fSmallestFraction= newFraction;
958 }
959 else
960 {
961 std::ostringstream message;
962 message << "Smallest Fraction not changed. " << G4endl
963 << " Proposed value was " << newFraction << G4endl
964 << " Value must be between 1.e-8 and 1.e-16";
965 G4Exception("G4MagInt_Driver::SetSmallestFraction()",
966 "GeomField1001", JustWarning, message);
967 }
968}

◆ SetVerboseLevel()

void G4MagInt_Driver::SetVerboseLevel ( G4int level)
overridevirtual

Setter and getter for verbosity.

Implements G4VIntegrationDriver.

◆ StreamInfo()

void G4MagInt_Driver::StreamInfo ( std::ostream & os) const
overridevirtual

Writes out to stream the parameters/state of the driver.

Implements G4VIntegrationDriver.

Definition at line 1016 of file G4MagIntegratorDriver.cc.

1017{
1018 os << "State of G4MagInt_Driver: " << std::endl;
1019 os << " Max number of Steps = " << fMaxNoSteps
1020 << " (base # = " << fMaxStepBase << " )" << std::endl;
1021 os << " Safety factor = " << safety << std::endl;
1022 os << " Power - shrink = " << pshrnk << std::endl;
1023 os << " Power - grow = " << pgrow << std::endl;
1024 os << " threshold (errcon) = " << errcon << std::endl;
1025
1026 os << " fMinimumStep = " << fMinimumStep << std::endl;
1027 os << " Smallest Fraction = " << fSmallestFraction << std::endl;
1028
1029 os << " No Integrat Vars = " << fNoIntegrationVariables << std::endl;
1030 os << " Min No Vars = " << fMinNoVars << std::endl;
1031 os << " Num-Vars = " << fNoVars << std::endl;
1032
1033 os << " verbose level = " << fVerboseLevel << std::endl;
1034 os << " Reintegrates = " << DoesReIntegrate() << std::endl;
1035}
G4bool DoesReIntegrate() const override

◆ WarnEndPointTooFar()

void G4MagInt_Driver::WarnEndPointTooFar ( G4double endPointDist,
G4double hStepSize,
G4double epsilonRelative,
G4int debugFlag )
protected

Definition at line 434 of file G4MagIntegratorDriver.cc.

438{
439 static G4ThreadLocal G4double maxRelError = 0.0;
440 G4bool isNewMax, prNewMax;
441
442 isNewMax = endPointDist > (1.0 + maxRelError) * h;
443 prNewMax = endPointDist > (1.0 + 1.05 * maxRelError) * h;
444 if( isNewMax ) { maxRelError= endPointDist / h - 1.0; }
445
446 if( (dbg != 0) && (h > G4GeometryTolerance::GetInstance()->GetSurfaceTolerance())
447 && ( (dbg>1) || prNewMax || (endPointDist >= h*(1.+eps) ) ) )
448 {
449 static G4ThreadLocal G4int noWarnings = 0;
450 std::ostringstream message;
451 if( (noWarnings++ < 10) || (dbg>2) )
452 {
453 message << "The integration produced an end-point which " << G4endl
454 << "is further from the start-point than the curve length."
455 << G4endl;
456 }
457 message << " Distance of endpoints = " << endPointDist
458 << ", curve length = " << h << G4endl
459 << " Difference (curveLen-endpDist)= " << (h - endPointDist)
460 << ", relative = " << (h-endPointDist) / h
461 << ", epsilon = " << eps;
462 G4Exception("G4MagInt_Driver::WarnEndPointTooFar()", "GeomField1001",
463 JustWarning, message);
464 }
465}
static G4GeometryTolerance * GetInstance()

Referenced by AccurateAdvance().

◆ WarnSmallStepSize()

void G4MagInt_Driver::WarnSmallStepSize ( G4double hnext,
G4double hstep,
G4double h,
G4double xDone,
G4int noSteps )
protected

Loggers, issuing warnings for undesirable situations.

Definition at line 385 of file G4MagIntegratorDriver.cc.

388{
389 static G4ThreadLocal G4int noWarningsIssued = 0;
390 const G4int maxNoWarnings = 10; // Number of verbose warnings
391 std::ostringstream message;
392 if( (noWarningsIssued < maxNoWarnings) || fVerboseLevel > 10 )
393 {
394 message << "The stepsize for the next iteration, " << hnext
395 << ", is too small - in Step number " << nstp << "." << G4endl
396 << "The minimum for the driver is " << Hmin() << G4endl
397 << "Requested integr. length was " << hstep << " ." << G4endl
398 << "The size of this sub-step was " << h << " ." << G4endl
399 << "The integrations has already gone " << xDone;
400 }
401 else
402 {
403 message << "Too small 'next' step " << hnext
404 << ", step-no: " << nstp << G4endl
405 << ", this sub-step: " << h
406 << ", req_tot_len: " << hstep
407 << ", done: " << xDone << ", min: " << Hmin();
408 }
409 G4Exception("G4MagInt_Driver::WarnSmallStepSize()", "GeomField1001",
410 JustWarning, message);
411 ++noWarningsIssued;
412}

Referenced by AccurateAdvance().

◆ WarnTooManyStep()

void G4MagInt_Driver::WarnTooManyStep ( G4double x1start,
G4double x2end,
G4double xCurrent )
protected

Definition at line 417 of file G4MagIntegratorDriver.cc.

420{
421 std::ostringstream message;
422 message << "The number of steps used in the Integration driver"
423 << " (Runge-Kutta) is too many." << G4endl
424 << "Integration of the interval was not completed !" << G4endl
425 << "Only a " << (xCurrent-x1start)*100/(x2end-x1start)
426 << " % fraction of it was done.";
427 G4Exception("G4MagInt_Driver::WarnTooManyStep()", "GeomField1001",
428 JustWarning, message);
429}

Referenced by AccurateAdvance().


The documentation for this class was generated from the following files: