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

#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 ()
 
void TestChordPrint (G4int noTrials, G4int lastStepTrial, G4double dChordStep, G4double fDeltaChord, G4double nextStepTrial)
 
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 ()
 
G4double GetLastStepEstimateUnc ()
 
void SetLastStepEstimateUnc (G4double stepEst)
 
void StreamDelegateInfo (std::ostream &os) const
 

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

Definition at line 44 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 )

Definition at line 49 of file G4MagIntegratorDriver.cc.

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

◆ ~G4MagInt_Driver()

G4MagInt_Driver::~G4MagInt_Driver ( )
override

Definition at line 85 of file G4MagIntegratorDriver.cc.

86{
87 if( fStatisticsVerboseLevel > 1 )
88 {
90 }
91}

◆ G4MagInt_Driver() [2/2]

G4MagInt_Driver::G4MagInt_Driver ( const G4MagInt_Driver & )
delete

Member Function Documentation

◆ AccurateAdvance()

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

Implements G4VIntegrationDriver.

Definition at line 96 of file G4MagIntegratorDriver.cc.

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

◆ AdvanceChordLimited()

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

Implements G4VIntegrationDriver.

◆ ComputeAndSetErrcon()

G4double G4MagInt_Driver::ComputeAndSetErrcon ( )
inline

◆ ComputeNewStepSize()

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

Implements G4VIntegrationDriver.

Definition at line 736 of file G4MagIntegratorDriver.cc.

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

Referenced by AccurateAdvance(), and QuickAdvance().

◆ ComputeNewStepSize_WithinLimits()

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

Definition at line 752 of file G4MagIntegratorDriver.cc.

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

Definition at line 707 of file G4MagIntegratorDriver.cc.

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

Referenced by ComputeNewStepSize().

◆ DoesReIntegrate()

G4bool G4MagInt_Driver::DoesReIntegrate ( ) const
inlineoverridevirtual

Implements G4VIntegrationDriver.

Definition at line 66 of file G4MagIntegratorDriver.hh.

66{ 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 981 of file G4MagIntegratorDriver.cc.

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

◆ GetDerivatives() [2/2]

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

Implements G4VIntegrationDriver.

◆ GetEquationOfMotion()

G4EquationOfMotion * G4MagInt_Driver::GetEquationOfMotion ( )
overridevirtual

Implements G4VIntegrationDriver.

Definition at line 990 of file G4MagIntegratorDriver.cc.

991{
992 return pIntStepper->GetEquationOfMotion();
993}
G4EquationOfMotion * GetEquationOfMotion()

◆ GetErrcon()

G4double G4MagInt_Driver::GetErrcon ( ) const
inline

Referenced by PrintInfo().

◆ GetHmin()

G4double G4MagInt_Driver::GetHmin ( ) const
inline

Referenced by PrintInfo().

◆ GetMaxNoSteps()

G4int G4MagInt_Driver::GetMaxNoSteps ( ) const
inline

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

Implements G4VIntegrationDriver.

Definition at line 1000 of file G4MagIntegratorDriver.cc.

1001{
1002 return pIntStepper;
1003}

◆ GetStepper() [2/2]

G4MagIntegratorStepper * G4MagInt_Driver::GetStepper ( )
overridevirtual

Implements G4VIntegrationDriver.

Definition at line 1005 of file G4MagIntegratorDriver.cc.

1006{
1007 return pIntStepper;
1008}

◆ 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

Implements G4VIntegrationDriver.

Definition at line 65 of file G4MagIntegratorDriver.hh.

65{}

◆ OneGoodStep()

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

Definition at line 471 of file G4MagIntegratorDriver.cc.

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

Referenced by AccurateAdvance().

◆ OnStartTracking()

void G4MagInt_Driver::OnStartTracking ( )
inlineoverridevirtual

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 872 of file G4MagIntegratorDriver.cc.

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

Definition at line 937 of file G4MagIntegratorDriver.cc.

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

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

Definition at line 784 of file G4MagIntegratorDriver.cc.

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

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 809 of file G4MagIntegratorDriver.cc.

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

Definition at line 595 of file G4MagIntegratorDriver.cc.

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

◆ QuickAdvance() [2/2]

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

Reimplemented from G4VIntegrationDriver.

Definition at line 614 of file G4MagIntegratorDriver.cc.

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

Referenced by AccurateAdvance().

◆ RenewStepperAndAdjust()

void G4MagInt_Driver::RenewStepperAndAdjust ( G4MagIntegratorStepper * pItsStepper)
overridevirtual

Reimplemented from G4VIntegrationDriver.

Definition at line 1010 of file G4MagIntegratorDriver.cc.

1012{
1013 pIntStepper = pItsStepper;
1015}
void ReSetParameters(G4double new_safety=0.9)

Referenced by G4MagInt_Driver().

◆ ReSetParameters()

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

Referenced by RenewStepperAndAdjust().

◆ SetEquationOfMotion()

void G4MagInt_Driver::SetEquationOfMotion ( G4EquationOfMotion * equation)
overridevirtual

Implements G4VIntegrationDriver.

Definition at line 995 of file G4MagIntegratorDriver.cc.

996{
997 pIntStepper->SetEquationOfMotion(equation);
998}
void SetEquationOfMotion(G4EquationOfMotion *newEquation)

◆ SetErrcon()

void G4MagInt_Driver::SetErrcon ( G4double valEc)
inline

◆ SetHmin()

void G4MagInt_Driver::SetHmin ( G4double newval)
inline

◆ 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

◆ SetSmallestFraction()

void G4MagInt_Driver::SetSmallestFraction ( G4double val)

Definition at line 954 of file G4MagIntegratorDriver.cc.

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

◆ SetVerboseLevel()

void G4MagInt_Driver::SetVerboseLevel ( G4int newLevel)
overridevirtual

Implements G4VIntegrationDriver.

◆ StreamInfo()

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

Implements G4VIntegrationDriver.

Definition at line 1017 of file G4MagIntegratorDriver.cc.

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

◆ WarnEndPointTooFar()

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

Definition at line 435 of file G4MagIntegratorDriver.cc.

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

Referenced by AccurateAdvance().

◆ WarnSmallStepSize()

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

Definition at line 386 of file G4MagIntegratorDriver.cc.

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

Referenced by AccurateAdvance().

◆ WarnTooManyStep()

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

Definition at line 418 of file G4MagIntegratorDriver.cc.

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

Referenced by AccurateAdvance().


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