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