Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MagIntegratorDriver.cc
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// G4MagInt_Driver implementation
27//
28// V.Grichine, 07.10.1996 - Created
29// W.Wander, 28.01.1998 - Added ability for low order integrators
30// J.Apostolakis, 08.11.2001 - Respect minimum step in AccurateAdvance
31// --------------------------------------------------------------------
32
33#include <iomanip>
34
35#include "globals.hh"
36#include "G4SystemOfUnits.hh"
39#include "G4FieldTrack.hh"
40
41#ifdef G4DEBUG_FIELD
42#include "G4DriverReporter.hh"
43#endif
44
45// ---------------------------------------------------------
46
47// Constructor
48//
50 G4MagIntegratorStepper* pStepper,
51 G4int numComponents,
52 G4int statisticsVerbose)
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}
80
81// ---------------------------------------------------------
82
83// Destructor
84//
86{
87 if( fStatisticsVerboseLevel > 1 )
88 {
90 }
91}
92
93// ---------------------------------------------------------
94
97 G4double hstep,
98 G4double eps,
99 G4double hinitial )
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, no_warnings = 0;
108 G4double x, hnext, hdid, h;
109
110#ifdef G4DEBUG_FIELD
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
119 G4double x1, x2;
120 G4bool succeeded = true, lastStepSucceeded;
121
122 G4double startCurveLength;
123
124 G4int noFullIntegr = 0, noSmallIntegr = 0;
125 static G4ThreadLocal G4int noGoodSteps = 0; // Bad = chord > curve-len
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 else
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
154 y_current.DumpToArray( ystart );
155
156 startCurveLength= y_current.GetCurveLength();
157 x1= startCurveLength;
158 x2= x1 + hstep;
159
160 if ( (hinitial > 0.0) && (hinitial < hstep)
161 && (hinitial > perMillion * hstep) )
162 {
163 h = hinitial;
164 }
165 else // Initial Step size "h" defaults to the full interval
166 {
167 h = hstep;
168 }
169
170 x = x1;
171
172 for ( i=0; i<nvar; ++i) { y[i] = ystart[i]; }
173
174 G4bool lastStep= false;
175 nstp = 1;
176
177 do
178 {
179 G4ThreeVector StartPos( y[0], y[1], y[2] );
180
181#ifdef G4DEBUG_FIELD
182 G4double xSubStepStart= x;
183 for (i=0; i<nvar; ++i) { ySubStepStart[i] = y[i]; }
184 yFldTrkStart.LoadFromArray(y, fNoIntegrationVariables);
185 yFldTrkStart.SetCurveLength(x);
186#endif
187
188 pIntStepper->RightHandSide( y, dydx );
189 ++fNoTotalSteps;
190
191 // Perform the Integration
192 //
193 if( h > fMinimumStep )
194 {
195 OneGoodStep(y,dydx,x,h,eps,hdid,hnext) ;
196 //--------------------------------------
197 lastStepSucceeded = (hdid == h);
198#ifdef G4DEBUG_FIELD
199 if (dbg>2)
200 {
201 // PrintStatus( ySubStepStart, xSubStepStart, y, x, h, nstp); // Only
202 G4DriverReporter::PrintStatus( ySubStepStart, xSubStepStart, y, x, h, nstp, nvar);
203 }
204#endif
205 }
206 else
207 {
208 G4FieldTrack yFldTrk( G4ThreeVector(0,0,0),
209 G4ThreeVector(0,0,0), 0., 0., 0., 0. );
210 G4double dchord_step, dyerr, dyerr_len; // What to do with these ?
211 yFldTrk.LoadFromArray(y, fNoIntegrationVariables);
212 yFldTrk.SetCurveLength( x );
213
214 QuickAdvance( yFldTrk, dydx, h, dchord_step, dyerr_len );
215 //-----------------------------------------------------
216
217 yFldTrk.DumpToArray(y);
218
219#ifdef G4FLD_STATS
220 ++fNoSmallSteps;
221 if ( dyerr_len > fDyerr_max ) { fDyerr_max = dyerr_len; }
222 fDyerrPos_smTot += dyerr_len;
223 fSumH_sm += h; // Length total for 'small' steps
224 if (nstp==1) { ++fNoInitialSmallSteps; }
225#endif
226#ifdef G4DEBUG_FIELD
227 if (dbg>1)
228 {
229 if(fNoSmallSteps<2) { PrintStatus(ySubStepStart, x1, y, x, h, -nstp); }
230 G4cout << "Another sub-min step, no " << fNoSmallSteps
231 << " of " << fNoTotalSteps << " this time " << nstp << G4endl;
232 PrintStatus( ySubStepStart, x1, y, x, h, nstp); // Only this
233 G4cout << " dyerr= " << dyerr_len << " relative = " << dyerr_len / h
234 << " epsilon= " << eps << " hstep= " << hstep
235 << " h= " << h << " hmin= " << fMinimumStep << G4endl;
236 }
237#endif
238 if( h == 0.0 )
239 {
240 G4Exception("G4MagInt_Driver::AccurateAdvance()",
241 "GeomField0003", FatalException,
242 "Integration Step became Zero!");
243 }
244 dyerr = dyerr_len / h;
245 hdid = h;
246 x += hdid;
247
248 // Compute suggested new step
249 hnext = ComputeNewStepSize( dyerr/eps, h);
250
251 // .. hnext= ComputeNewStepSize_WithinLimits( dyerr/eps, h);
252 lastStepSucceeded = (dyerr<= eps);
253 }
254
255 if (lastStepSucceeded) { ++noFullIntegr; }
256 else { ++noSmallIntegr; }
257
258 G4ThreeVector EndPos( y[0], y[1], y[2] );
259
260#ifdef G4DEBUG_FIELD
261 if( (dbg>0) && (dbg<=2) && (nstp>nStpPr))
262 {
263 if( nstp==nStpPr ) { G4cout << "***** Many steps ****" << G4endl; }
264 G4cout << "MagIntDrv: " ;
265 G4cout << "hdid=" << std::setw(12) << hdid << " "
266 << "hnext=" << std::setw(12) << hnext << " "
267 << "hstep=" << std::setw(12) << hstep << " (requested) "
268 << G4endl;
269 PrintStatus( ystart, x1, y, x, h, (nstp==nStpPr) ? -nstp: nstp);
270 }
271#endif
272
273 // Check the endpoint
274 G4double endPointDist= (EndPos-StartPos).mag();
275 if ( endPointDist >= hdid*(1.+perMillion) )
276 {
277 ++fNoBadSteps;
278
279 // Issue a warning only for gross differences -
280 // we understand how small difference occur.
281 if ( endPointDist >= hdid*(1.+perThousand) )
282 {
283#ifdef G4DEBUG_FIELD
284 if (dbg)
285 {
286 WarnEndPointTooFar ( endPointDist, hdid, eps, dbg );
287 G4cerr << " Total steps: bad " << fNoBadSteps
288 << " good " << noGoodSteps << " current h= " << hdid
289 << G4endl;
290 PrintStatus( ystart, x1, y, x, hstep, no_warnings?nstp:-nstp);
291 }
292#endif
293 ++no_warnings;
294 }
295 }
296 else
297 {
298 ++noGoodSteps;
299 }
300// #endif
301
302 // Avoid numerous small last steps
303 if( (h < eps * hstep) || (h < fSmallestFraction * startCurveLength) )
304 {
305 // No more integration -- the next step will not happen
306 lastStep = true;
307 }
308 else
309 {
310 // Check the proposed next stepsize
311 if(std::fabs(hnext) <= Hmin())
312 {
313#ifdef G4DEBUG_FIELD
314 // If simply a very small interval is being integrated, do not warn
315 if( (x < x2 * (1-eps) ) && // The last step can be small: OK
316 (std::fabs(hstep) > Hmin()) ) // and if we are asked, it's OK
317 {
318 if(dbg>0)
319 {
320 WarnSmallStepSize( hnext, hstep, h, x-x1, nstp );
321 PrintStatus( ystart, x1, y, x, hstep, no_warnings?nstp:-nstp);
322 }
323 ++no_warnings;
324 }
325#endif
326 // Make sure that the next step is at least Hmin.
327 h = Hmin();
328 }
329 else
330 {
331 h = hnext;
332 }
333
334 // Ensure that the next step does not overshoot
335 if ( x+h > x2 )
336 { // When stepsize overshoots, decrease it!
337 h = x2 - x ; // Must cope with difficult rounding-error
338 } // issues if hstep << x2
339
340 if ( h == 0.0 )
341 {
342 // Cannot progress - accept this as last step - by default
343 lastStep = true;
344#ifdef G4DEBUG_FIELD
345 if (dbg>2)
346 {
347 int prec= G4cout.precision(12);
348 G4cout << "Warning: G4MagIntegratorDriver::AccurateAdvance"
349 << G4endl
350 << " Integration step 'h' became "
351 << h << " due to roundoff. " << G4endl
352 << " Calculated as difference of x2= "<< x2 << " and x=" << x
353 << " Forcing termination of advance." << G4endl;
354 G4cout.precision(prec);
355 }
356#endif
357 }
358 }
359 } while ( ((++nstp)<=fMaxNoSteps) && (x < x2) && (!lastStep) );
360 // Loop checking, 07.10.2016, J. Apostolakis
361
362 // Have we reached the end ?
363 // --> a better test might be x-x2 > an_epsilon
364
365 succeeded = (x>=x2); // If it was a "forced" last step
366
367 for (i=0; i<nvar; ++i) { yEnd[i] = y[i]; }
368
369 // Put back the values.
370 y_current.LoadFromArray( yEnd, fNoIntegrationVariables );
371 y_current.SetCurveLength( x );
372
373 if(nstp > fMaxNoSteps)
374 {
375 ++no_warnings;
376 succeeded = false;
377#ifdef G4DEBUG_FIELD
378 if (dbg)
379 {
380 WarnTooManyStep( x1, x2, x ); // Issue WARNING
381 PrintStatus( yEnd, x1, y, x, hstep, -nstp);
382 }
383#endif
384 }
385
386#ifdef G4DEBUG_FIELD
387 if( dbg && no_warnings )
388 {
389 G4cerr << "G4MagIntegratorDriver exit status: no-steps " << nstp << G4endl;
390 PrintStatus( yEnd, x1, y, x, hstep, nstp);
391 }
392#endif
393
394 return succeeded;
395} // end of AccurateAdvance ...........................
396
397// ---------------------------------------------------------
398
399void
401 G4double h, G4double xDone,
402 G4int nstp)
403{
404 static G4ThreadLocal G4int noWarningsIssued = 0;
405 const G4int maxNoWarnings = 10; // Number of verbose warnings
406 std::ostringstream message;
407 if( (noWarningsIssued < maxNoWarnings) || fVerboseLevel > 10 )
408 {
409 message << "The stepsize for the next iteration, " << hnext
410 << ", is too small - in Step number " << nstp << "." << G4endl
411 << "The minimum for the driver is " << Hmin() << G4endl
412 << "Requested integr. length was " << hstep << " ." << G4endl
413 << "The size of this sub-step was " << h << " ." << G4endl
414 << "The integrations has already gone " << xDone;
415 }
416 else
417 {
418 message << "Too small 'next' step " << hnext
419 << ", step-no: " << nstp << G4endl
420 << ", this sub-step: " << h
421 << ", req_tot_len: " << hstep
422 << ", done: " << xDone << ", min: " << Hmin();
423 }
424 G4Exception("G4MagInt_Driver::WarnSmallStepSize()", "GeomField1001",
425 JustWarning, message);
426 ++noWarningsIssued;
427}
428
429// ---------------------------------------------------------
430
431void
433 G4double x2end,
434 G4double xCurrent )
435{
436 std::ostringstream message;
437 message << "The number of steps used in the Integration driver"
438 << " (Runge-Kutta) is too many." << G4endl
439 << "Integration of the interval was not completed !" << G4endl
440 << "Only a " << (xCurrent-x1start)*100/(x2end-x1start)
441 << " % fraction of it was done.";
442 G4Exception("G4MagInt_Driver::WarnTooManyStep()", "GeomField1001",
443 JustWarning, message);
444}
445
446// ---------------------------------------------------------
447
448void
450 G4double h ,
451 G4double eps,
452 G4int dbg)
453{
454 static G4ThreadLocal G4double maxRelError = 0.0;
455 G4bool isNewMax, prNewMax;
456
457 isNewMax = endPointDist > (1.0 + maxRelError) * h;
458 prNewMax = endPointDist > (1.0 + 1.05 * maxRelError) * h;
459 if( isNewMax ) { maxRelError= endPointDist / h - 1.0; }
460
462 && ( (dbg>1) || prNewMax || (endPointDist >= h*(1.+eps) ) ) )
463 {
464 static G4ThreadLocal G4int noWarnings = 0;
465 std::ostringstream message;
466 if( (noWarnings++ < 10) || (dbg>2) )
467 {
468 message << "The integration produced an end-point which " << G4endl
469 << "is further from the start-point than the curve length."
470 << G4endl;
471 }
472 message << " Distance of endpoints = " << endPointDist
473 << ", curve length = " << h << G4endl
474 << " Difference (curveLen-endpDist)= " << (h - endPointDist)
475 << ", relative = " << (h-endPointDist) / h
476 << ", epsilon = " << eps;
477 G4Exception("G4MagInt_Driver::WarnEndPointTooFar()", "GeomField1001",
478 JustWarning, message);
479 }
480}
481
482// ---------------------------------------------------------
483
484void
486 const G4double dydx[],
487 G4double& x, // InOut
488 G4double htry,
489 G4double eps_rel_max,
490 G4double& hdid, // Out
491 G4double& hnext ) // Out
492
493// Driver for one Runge-Kutta Step with monitoring of local truncation error
494// to ensure accuracy and adjust stepsize. Input are dependent variable
495// array y[0,...,5] and its derivative dydx[0,...,5] at the
496// starting value of the independent variable x . Also input are stepsize
497// to be attempted htry, and the required accuracy eps. On output y and x
498// are replaced by their new values, hdid is the stepsize that was actually
499// accomplished, and hnext is the estimated next stepsize.
500// This is similar to the function rkqs from the book:
501// Numerical Recipes in C: The Art of Scientific Computing (NRC), Second
502// Edition, by William H. Press, Saul A. Teukolsky, William T.
503// Vetterling, and Brian P. Flannery (Cambridge University Press 1992),
504// 16.2 Adaptive StepSize Control for Runge-Kutta, p. 719
505
506{
507 G4double errmax_sq;
508 G4double h, htemp, xnew ;
509
511
512 h = htry ; // Set stepsize to the initial trial value
513
514 G4double inv_eps_vel_sq = 1.0 / (eps_rel_max*eps_rel_max);
515
516 G4double errpos_sq = 0.0; // square of displacement error
517 G4double errvel_sq = 0.0; // square of momentum vector difference
518 G4double errspin_sq = 0.0; // square of spin vector difference
519
520 static G4ThreadLocal G4int tot_no_trials=0;
521 const G4int max_trials=100;
522
523 G4ThreeVector Spin(y[9],y[10],y[11]);
524 G4double spin_mag2 = Spin.mag2();
525 G4bool hasSpin = (spin_mag2 > 0.0);
526
527 for (G4int iter=0; iter<max_trials; ++iter)
528 {
529 ++tot_no_trials;
530 pIntStepper-> Stepper(y,dydx,h,ytemp,yerr);
531 // *******
532 G4double eps_pos = eps_rel_max * std::max(h, fMinimumStep);
533 G4double inv_eps_pos_sq = 1.0 / (eps_pos*eps_pos);
534
535 // Evaluate accuracy
536 //
537 errpos_sq = sqr(yerr[0]) + sqr(yerr[1]) + sqr(yerr[2]) ;
538 errpos_sq *= inv_eps_pos_sq; // Scale relative to required tolerance
539
540 // Accuracy for momentum
541 G4double magvel_sq= sqr(y[3]) + sqr(y[4]) + sqr(y[5]) ;
542 G4double sumerr_sq = sqr(yerr[3]) + sqr(yerr[4]) + sqr(yerr[5]) ;
543 if( magvel_sq > 0.0 )
544 {
545 errvel_sq = sumerr_sq / magvel_sq;
546 }
547 else
548 {
549 std::ostringstream message;
550 message << "Found case of zero momentum." << G4endl
551 << "- iteration= " << iter << "; h= " << h;
552 G4Exception("G4MagInt_Driver::OneGoodStep()",
553 "GeomField1001", JustWarning, message);
554 errvel_sq = sumerr_sq;
555 }
556 errvel_sq *= inv_eps_vel_sq;
557 errmax_sq = std::max( errpos_sq, errvel_sq ); // Square of maximum error
558
559 if( hasSpin )
560 {
561 // Accuracy for spin
562 errspin_sq = ( sqr(yerr[9]) + sqr(yerr[10]) + sqr(yerr[11]) )
563 / spin_mag2; // ( sqr(y[9]) + sqr(y[10]) + sqr(y[11]) );
564 errspin_sq *= inv_eps_vel_sq;
565 errmax_sq = std::max( errmax_sq, errspin_sq );
566 }
567
568 if ( errmax_sq <= 1.0 ) { break; } // Step succeeded.
569
570 // Step failed; compute the size of retrial Step.
571 htemp = GetSafety() * h * std::pow( errmax_sq, 0.5*GetPshrnk() );
572
573 if (htemp >= 0.1*h) { h = htemp; } // Truncation error too large,
574 else { h = 0.1*h; } // reduce stepsize, but no more
575 // than a factor of 10
576 xnew = x + h;
577 if(xnew == x)
578 {
579 std::ostringstream message;
580 message << "Stepsize underflow in Stepper !" << G4endl
581 << "- Step's start x=" << x << " and end x= " << xnew
582 << " are equal !! " << G4endl
583 << " Due to step-size= " << h
584 << ". Note that input step was " << htry;
585 G4Exception("G4MagInt_Driver::OneGoodStep()",
586 "GeomField1001", JustWarning, message);
587 break;
588 }
589 }
590
591 // Compute size of next Step
592 if (errmax_sq > errcon*errcon)
593 {
594 hnext = GetSafety()*h*std::pow(errmax_sq, 0.5*GetPgrow());
595 }
596 else
597 {
598 hnext = max_stepping_increase*h ; // No more than a factor of 5 increase
599 }
600 x += (hdid = h);
601
602 for(G4int k=0; k<fNoIntegrationVariables; ++k) { y[k] = ytemp[k]; }
603
604 return;
605}
606
607//----------------------------------------------------------------------
608
609// QuickAdvance just tries one Step - it does not ensure accuracy
610//
612 const G4double dydx[],
613 G4double hstep, // In
614 G4double& dchord_step,
615 G4double& dyerr_pos_sq,
616 G4double& dyerr_mom_rel_sq )
617{
618 G4Exception("G4MagInt_Driver::QuickAdvance()", "GeomField0001",
619 FatalException, "Not yet implemented.");
620
621 // Use the parameters of this method, to please compiler
622 //
623 dchord_step = dyerr_pos_sq = hstep * hstep * dydx[0];
624 dyerr_mom_rel_sq = y_posvel.GetPosition().mag2();
625 return true;
626}
627
628//----------------------------------------------------------------------
629
631 const G4double dydx[],
632 G4double hstep, // In
633 G4double& dchord_step,
634 G4double& dyerr )
635{
636 G4double dyerr_pos_sq, dyerr_mom_rel_sq;
639 G4double s_start;
640 G4double dyerr_mom_sq, vel_mag_sq, inv_vel_mag_sq;
641
642 static G4ThreadLocal G4int no_call = 0;
643 ++no_call;
644
645 // Move data into array
646 y_posvel.DumpToArray( yarrin ); // yarrin <== y_posvel
647 s_start = y_posvel.GetCurveLength();
648
649 // Do an Integration Step
650 pIntStepper-> Stepper(yarrin, dydx, hstep, yarrout, yerr_vec) ;
651
652 // Estimate curve-chord distance
653 dchord_step= pIntStepper-> DistChord();
654
655 // Put back the values. yarrout ==> y_posvel
656 y_posvel.LoadFromArray( yarrout, fNoIntegrationVariables );
657 y_posvel.SetCurveLength( s_start + hstep );
658
659#ifdef G4DEBUG_FIELD
660 if(fVerboseLevel>2)
661 {
662 G4cout << "G4MagIntDrv: Quick Advance" << G4endl;
663 PrintStatus( yarrin, s_start, yarrout, s_start+hstep, hstep, 1);
664 }
665#endif
666
667 // A single measure of the error
668 // TO-DO : account for energy, spin, ... ?
669 vel_mag_sq = ( sqr(yarrout[3])+sqr(yarrout[4])+sqr(yarrout[5]) );
670 inv_vel_mag_sq = 1.0 / vel_mag_sq;
671 dyerr_pos_sq = ( sqr(yerr_vec[0])+sqr(yerr_vec[1])+sqr(yerr_vec[2]));
672 dyerr_mom_sq = ( sqr(yerr_vec[3])+sqr(yerr_vec[4])+sqr(yerr_vec[5]));
673 dyerr_mom_rel_sq = dyerr_mom_sq * inv_vel_mag_sq;
674
675 // Calculate also the change in the momentum squared also ???
676 // G4double veloc_square = y_posvel.GetVelocity().mag2();
677 // ...
678
679#ifdef RETURN_A_NEW_STEP_LENGTH
680 // The following step cannot be done here because "eps" is not known.
681 dyerr_len = std::sqrt( dyerr_len_sq );
682 dyerr_len_sq /= eps ;
683
684 // Look at the velocity deviation ?
685 // sqr(yerr_vec[3])+sqr(yerr_vec[4])+sqr(yerr_vec[5]));
686
687 // Set suggested new step
688 hstep = ComputeNewStepSize( dyerr_len, hstep);
689#endif
690
691 if( dyerr_pos_sq > ( dyerr_mom_rel_sq * sqr(hstep) ) )
692 {
693 dyerr = std::sqrt(dyerr_pos_sq);
694 }
695 else
696 {
697 // Scale it to the current step size - for now
698 dyerr = std::sqrt(dyerr_mom_rel_sq) * hstep;
699 }
700
701 return true;
702}
703
704// --------------------------------------------------------------------------
705
706#ifdef QUICK_ADV_ARRAY_IN_AND_OUT
708 const G4double dydx[],
709 G4double hstep, // In
710 G4double yarrout[],
711 G4double& dchord_step,
712 G4double& dyerr ) // In length
713{
714 G4Exception("G4MagInt_Driver::QuickAdvance()", "GeomField0001",
715 FatalException, "Not yet implemented.");
716 dyerr = dchord_step = hstep * yarrin[0] * dydx[0];
717 yarrout[0]= yarrin[0];
718}
719#endif
720
721// --------------------------------------------------------------------------
722
723// This method computes new step sizes - but does not limit changes to
724// within certain factors
725//
727ComputeNewStepSize_WithoutReductionLimit(G4double errMaxNorm, // max error (normalised)
728 G4double hstepCurrent) // current step size
729{
730 G4double hnew;
731
732 // Compute size of next Step for a failed step
733 if(errMaxNorm > 1.0 )
734 {
735 // Step failed; compute the size of retrial Step.
736 hnew = GetSafety()*hstepCurrent*std::pow(errMaxNorm,GetPshrnk()) ;
737 }
738 else if(errMaxNorm > 0.0 )
739 {
740 // Compute size of next Step for a successful step
741 hnew = GetSafety()*hstepCurrent*std::pow(errMaxNorm,GetPgrow()) ;
742 }
743 else
744 {
745 // if error estimate is zero (possible) or negative (dubious)
746 hnew = max_stepping_increase * hstepCurrent;
747 }
748
749 return hnew;
750}
751
752// ---------------------------------------------------------------------------
753
756 G4double errMaxNorm, // max error (normalised)
757 G4double hstepCurrent) // current step size
758{
759 // Legacy behaviour:
760 return ComputeNewStepSize_WithoutReductionLimit( errMaxNorm, hstepCurrent );
761 // 'Improved' behaviour - at least more consistent with other step estimates:
762 // return ComputeNewStepSize_WithinLimits( errMaxNorm, hstepCurrent );
763}
764
765// This method computes new step sizes limiting changes within certain factors
766//
767// It shares its logic with AccurateAdvance.
768// They are kept separate currently for optimisation.
769//
772 G4double errMaxNorm, // max error (normalised)
773 G4double hstepCurrent) // current step size
774{
775 G4double hnew;
776
777 // Compute size of next Step for a failed step
778 if (errMaxNorm > 1.0 )
779 {
780 // Step failed; compute the size of retrial Step.
781 hnew = GetSafety()*hstepCurrent*std::pow(errMaxNorm,GetPshrnk()) ;
782
783 if (hnew < max_stepping_decrease*hstepCurrent)
784 {
785 hnew = max_stepping_decrease*hstepCurrent ;
786 // reduce stepsize, but no more
787 // than this factor (value= 1/10)
788 }
789 }
790 else
791 {
792 // Compute size of next Step for a successful step
793 if (errMaxNorm > errcon)
794 { hnew = GetSafety()*hstepCurrent*std::pow(errMaxNorm,GetPgrow()); }
795 else // No more than a factor of 5 increase
796 { hnew = max_stepping_increase * hstepCurrent; }
797 }
798 return hnew;
799}
800
801// ---------------------------------------------------------------------------
802
804 G4double xstart,
805 const G4double* CurrentArr,
806 G4double xcurrent,
807 G4double requestStep,
808 G4int subStepNo )
809 // Potentially add as arguments:
810 // <dydx> - as Initial Force
811 // stepTaken(hdid) - last step taken
812 // nextStep (hnext) - proposal for size
813{
814 G4FieldTrack StartFT(G4ThreeVector(0,0,0),
815 G4ThreeVector(0,0,0), 0., 0., 0., 0. );
816 G4FieldTrack CurrentFT (StartFT);
817
818 StartFT.LoadFromArray( StartArr, fNoIntegrationVariables);
819 StartFT.SetCurveLength( xstart);
820 CurrentFT.LoadFromArray( CurrentArr, fNoIntegrationVariables);
821 CurrentFT.SetCurveLength( xcurrent );
822
823 PrintStatus(StartFT, CurrentFT, requestStep, subStepNo );
824}
825
826// ---------------------------------------------------------------------------
827
829 const G4FieldTrack& CurrentFT,
830 G4double requestStep,
831 G4int subStepNo)
832{
833 G4int verboseLevel= fVerboseLevel;
834 const G4int noPrecision = 5;
835 G4int oldPrec= G4cout.precision(noPrecision);
836 // G4cout.setf(ios_base::fixed,ios_base::floatfield);
837
838 const G4ThreeVector StartPosition= StartFT.GetPosition();
839 const G4ThreeVector StartUnitVelocity= StartFT.GetMomentumDir();
840 const G4ThreeVector CurrentPosition= CurrentFT.GetPosition();
841 const G4ThreeVector CurrentUnitVelocity= CurrentFT.GetMomentumDir();
842
843 G4double DotStartCurrentVeloc= StartUnitVelocity.dot(CurrentUnitVelocity);
844
845 G4double step_len= CurrentFT.GetCurveLength() - StartFT.GetCurveLength();
846 G4double subStepSize = step_len;
847
848 if( (subStepNo <= 1) || (verboseLevel > 3) )
849 {
850 subStepNo = - subStepNo; // To allow printing banner
851
852 G4cout << std::setw( 6) << " " << std::setw( 25)
853 << " G4MagInt_Driver: Current Position and Direction" << " "
854 << G4endl;
855 G4cout << std::setw( 5) << "Step#" << " "
856 << std::setw( 7) << "s-curve" << " "
857 << std::setw( 9) << "X(mm)" << " "
858 << std::setw( 9) << "Y(mm)" << " "
859 << std::setw( 9) << "Z(mm)" << " "
860 << std::setw( 8) << " N_x " << " "
861 << std::setw( 8) << " N_y " << " "
862 << std::setw( 8) << " N_z " << " "
863 << std::setw( 8) << " N^2-1 " << " "
864 << std::setw(10) << " N(0).N " << " "
865 << std::setw( 7) << "KinEner " << " "
866 << std::setw(12) << "Track-l" << " " // Add the Sub-step ??
867 << std::setw(12) << "Step-len" << " "
868 << std::setw(12) << "Step-len" << " "
869 << std::setw( 9) << "ReqStep" << " "
870 << G4endl;
871 }
872
873 if( (subStepNo <= 0) )
874 {
875 PrintStat_Aux( StartFT, requestStep, 0.,
876 0, 0.0, 1.0);
877 }
878
879 if( verboseLevel <= 3 )
880 {
881 G4cout.precision(noPrecision);
882 PrintStat_Aux( CurrentFT, requestStep, step_len,
883 subStepNo, subStepSize, DotStartCurrentVeloc );
884 }
885
886 G4cout.precision(oldPrec);
887}
888
889// ---------------------------------------------------------------------------
890
892 G4double requestStep,
893 G4double step_len,
894 G4int subStepNo,
895 G4double subStepSize,
896 G4double dotVeloc_StartCurr)
897{
898 const G4ThreeVector Position = aFieldTrack.GetPosition();
899 const G4ThreeVector UnitVelocity = aFieldTrack.GetMomentumDir();
900
901 if( subStepNo >= 0)
902 {
903 G4cout << std::setw( 5) << subStepNo << " ";
904 }
905 else
906 {
907 G4cout << std::setw( 5) << "Start" << " ";
908 }
909 G4double curveLen= aFieldTrack.GetCurveLength();
910 G4cout << std::setw( 7) << curveLen;
911 G4cout << std::setw( 9) << Position.x() << " "
912 << std::setw( 9) << Position.y() << " "
913 << std::setw( 9) << Position.z() << " "
914 << std::setw( 8) << UnitVelocity.x() << " "
915 << std::setw( 8) << UnitVelocity.y() << " "
916 << std::setw( 8) << UnitVelocity.z() << " ";
917 G4int oldprec= G4cout.precision(3);
918 G4cout << std::setw( 8) << UnitVelocity.mag2()-1.0 << " ";
919 G4cout.precision(6);
920 G4cout << std::setw(10) << dotVeloc_StartCurr << " ";
921 G4cout.precision(oldprec);
922 G4cout << std::setw( 7) << aFieldTrack.GetKineticEnergy();
923 G4cout << std::setw(12) << step_len << " ";
924
925 static G4ThreadLocal G4double oldCurveLength = 0.0;
926 static G4ThreadLocal G4double oldSubStepLength = 0.0;
927 static G4ThreadLocal G4int oldSubStepNo = -1;
928
929 G4double subStep_len = 0.0;
930 if( curveLen > oldCurveLength )
931 {
932 subStep_len= curveLen - oldCurveLength;
933 }
934 else if (subStepNo == oldSubStepNo)
935 {
936 subStep_len= oldSubStepLength;
937 }
938 oldCurveLength= curveLen;
939 oldSubStepLength= subStep_len;
940
941 G4cout << std::setw(12) << subStep_len << " ";
942 G4cout << std::setw(12) << subStepSize << " ";
943 if( requestStep != -1.0 )
944 {
945 G4cout << std::setw( 9) << requestStep << " ";
946 }
947 else
948 {
949 G4cout << std::setw( 9) << " InitialStep " << " ";
950 }
951 G4cout << G4endl;
952}
953
954// ---------------------------------------------------------------------------
955
957{
958 G4int noPrecBig = 6;
959 G4int oldPrec = G4cout.precision(noPrecBig);
960
961 G4cout << "G4MagInt_Driver Statistics of steps undertaken. " << G4endl;
962 G4cout << "G4MagInt_Driver: Number of Steps: "
963 << " Total= " << fNoTotalSteps
964 << " Bad= " << fNoBadSteps
965 << " Small= " << fNoSmallSteps
966 << " Non-initial small= " << (fNoSmallSteps-fNoInitialSmallSteps)
967 << G4endl;
968 G4cout.precision(oldPrec);
969}
970
971// ---------------------------------------------------------------------------
972
974{
975 if( (newFraction > 1.e-16) && (newFraction < 1e-8) )
976 {
977 fSmallestFraction= newFraction;
978 }
979 else
980 {
981 std::ostringstream message;
982 message << "Smallest Fraction not changed. " << G4endl
983 << " Proposed value was " << newFraction << G4endl
984 << " Value must be between 1.e-8 and 1.e-16";
985 G4Exception("G4MagInt_Driver::SetSmallestFraction()",
986 "GeomField1001", JustWarning, message);
987 }
988}
989
991GetDerivatives(const G4FieldTrack& y_curr, G4double* dydx) const
992{
994 y_curr.DumpToArray(ytemp);
995 pIntStepper->RightHandSide(ytemp, dydx);
996 // Avoid virtual call for GetStepper
997 // Was: GetStepper()->ComputeRightHandSide(ytemp, dydx);
998}
999
1001 G4double dydx[],
1002 G4double field[]) const
1003{
1005 track.DumpToArray(ytemp);
1006 pIntStepper->RightHandSide(ytemp, dydx, field);
1007}
1008
1010{
1011 return pIntStepper->GetEquationOfMotion();
1012}
1013
1015{
1016 pIntStepper->SetEquationOfMotion(equation);
1017}
1018
1020{
1021 return pIntStepper;
1022}
1023
1025{
1026 return pIntStepper;
1027}
1028
1031{
1032 pIntStepper = pItsStepper;
1034}
1035
1036void G4MagInt_Driver::StreamInfo( std::ostream& os ) const
1037{
1038 os << "State of G4MagInt_Driver: " << std::endl;
1039 os << " Max number of Steps = " << fMaxNoSteps
1040 << " (base # = " << fMaxStepBase << " )" << std::endl;
1041 os << " Safety factor = " << safety << std::endl;
1042 os << " Power - shrink = " << pshrnk << std::endl;
1043 os << " Power - grow = " << pgrow << std::endl;
1044 os << " threshold (errcon) = " << errcon << std::endl;
1045
1046 os << " fMinimumStep = " << fMinimumStep << std::endl;
1047 os << " Smallest Fraction = " << fSmallestFraction << std::endl;
1048
1049 os << " No Integrat Vars = " << fNoIntegrationVariables << std::endl;
1050 os << " Min No Vars = " << fMinNoVars << std::endl;
1051 os << " Num-Vars = " << fNoVars << std::endl;
1052
1053 os << " verbose level = " << fVerboseLevel << std::endl;
1054 os << " Reintegrates = " << DoesReIntegrate() << std::endl;
1055}
1056
1057void PrintInfo( const G4MagInt_Driver & magDrv, std::ostream& os )
1058{
1059 os << "State of G4MagInt_Driver: " << std::endl;
1060 os << " Max number of Steps = " << magDrv.GetMaxNoSteps();
1061 // << " (base # = " << magDrv.fMaxStepBase << " )" << std::endl;
1062 os << " Safety factor = " << magDrv.GetSafety() << std::endl;
1063 os << " Power - shrink = " << magDrv.GetPshrnk() << std::endl;
1064 os << " Power - grow = " << magDrv.GetPgrow() << std::endl;
1065 os << " threshold (errcon) = " << magDrv.GetErrcon() << std::endl;
1066
1067 os << " fMinimumStep = " << magDrv.GetHmin() << std::endl;
1068 os << " Smallest Fraction = " << magDrv.GetSmallestFraction() << std::endl;
1069
1070 /*****
1071 os << " No Integrat Vars = " << magDrv.GetNoIntegrationVariables << std::endl;
1072 os << " Min No Vars = " << magDrv.GetMinNoVars << std::endl;
1073 os << " Num-Vars = " << magDrv.GetNoVars << std::endl;
1074 *****/
1075 os << " verbose level = " << magDrv.GetVerboseLevel() << std::endl;
1076 os << " Reintegrates = " << magDrv.DoesReIntegrate() << std::endl;
1077}
const G4int noPrecision
@ JustWarning
@ FatalException
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
void PrintInfo(const G4MagInt_Driver &magDrv, std::ostream &os)
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
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double mag2() const
double y() const
double dot(const Hep3Vector &) const
static void PrintStatus(const G4double *StartArr, G4double xstart, const G4double *CurrentArr, G4double xcurrent, G4double requestStep, unsigned int subStepNo, unsigned int noIntegrationVariables)
const G4ThreeVector & GetMomentumDir() const
G4double GetCurveLength() const
void SetCurveLength(G4double nCurve_s)
G4double GetKineticEnergy() const
G4ThreeVector GetPosition() const
void DumpToArray(G4double valArr[ncompSVEC]) const
void LoadFromArray(const G4double valArr[ncompSVEC], G4int noVarsIntegrated)
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4double GetPshrnk() const
G4double ComputeNewStepSize_WithinLimits(G4double errMaxNorm, G4double hstepCurrent)
virtual void GetDerivatives(const G4FieldTrack &y_curr, G4double dydx[]) const override
G4MagInt_Driver(G4double hminimum, G4MagIntegratorStepper *pItsStepper, G4int numberOfComponents=6, G4int statisticsVerbosity=0)
void PrintStatus(const G4double *StartArr, G4double xstart, const G4double *CurrentArr, G4double xcurrent, G4double requestStep, G4int subStepNo)
virtual ~G4MagInt_Driver() override
virtual void SetEquationOfMotion(G4EquationOfMotion *equation) override
virtual void RenewStepperAndAdjust(G4MagIntegratorStepper *pItsStepper) override
G4double GetErrcon() const
void SetSmallestFraction(G4double val)
virtual const G4MagIntegratorStepper * GetStepper() const override
virtual G4bool QuickAdvance(G4FieldTrack &y_val, const G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr) override
G4double GetSafety() const
G4double Hmin() const
virtual G4bool AccurateAdvance(G4FieldTrack &y_current, G4double hstep, G4double eps, G4double hinitial=0.0) override
void WarnSmallStepSize(G4double hnext, G4double hstep, G4double h, G4double xDone, G4int noSteps)
virtual G4int GetVerboseLevel() const override
void WarnEndPointTooFar(G4double endPointDist, G4double hStepSize, G4double epsilonRelative, G4int debugFlag)
G4int GetMaxNoSteps() const
void OneGoodStep(G4double ystart[], const G4double dydx[], G4double &x, G4double htry, G4double eps, G4double &hdid, G4double &hnext)
G4double GetHmin() const
G4double ComputeNewStepSize_WithoutReductionLimit(G4double errMaxNorm, G4double hstepCurrent)
virtual G4EquationOfMotion * GetEquationOfMotion() override
void PrintStat_Aux(const G4FieldTrack &aFieldTrack, G4double requestStep, G4double actualStep, G4int subStepNo, G4double subStepSize, G4double dotVelocities)
G4double GetSmallestFraction() const
G4double GetPgrow() const
virtual G4bool DoesReIntegrate() const override
virtual G4double ComputeNewStepSize(G4double errMaxNorm, G4double hstepCurrent) override
void WarnTooManyStep(G4double x1start, G4double x2end, G4double xCurrent)
void StreamInfo(std::ostream &os) const override
void ReSetParameters(G4double new_safety=0.9)
G4EquationOfMotion * GetEquationOfMotion()
void RightHandSide(const G4double y[], G4double dydx[]) const
void SetEquationOfMotion(G4EquationOfMotion *newEquation)
virtual G4int IntegratorOrder() const =0
static constexpr G4double max_stepping_decrease
static constexpr G4double max_stepping_increase
T sqr(const T &x)
Definition: templates.hh:128
#define G4ThreadLocal
Definition: tls.hh:77