Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TDormandPrince45.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4TDormandPrince45
27//
28// Class desription:
29//
30// An implementation of the 5th order embedded RK method from the paper:
31// J. R. Dormand and P. J. Prince, "A family of embedded Runge-Kutta formulae"
32// Journal of computational and applied Math., vol.6, no.1, pp.19-26, 1980.
33//
34// DormandPrince7 - 5(4) embedded RK method
35//
36
37// Created: Somnath Banerjee, Google Summer of Code 2015, 25 May 2015
38// Supervision: John Apostolakis, CERN
39// --------------------------------------------------------------------
40#ifndef G4TDORMAND_PRINCE_45_HH
41#define G4TDORMAND_PRINCE_45_HH
42
44#include "G4FieldUtils.hh"
45#include "G4LineSection.hh"
46
47#include <cstring>
48#include <cassert>
49
50template <class T_Equation, unsigned int N = 6 >
52{
53 public:
54
55 G4TDormandPrince45(T_Equation* equation );
56 G4TDormandPrince45(T_Equation* equation, G4int numVar ); // must have numVar == N
57
58 inline void StepWithError(const G4double yInput[],
59 const G4double dydx[],
60 G4double hstep,
61 G4double yOutput[],
62 G4double yError[] ) ;
63
64 void Stepper(const G4double yInput[],
65 const G4double dydx[],
66 G4double hstep,
67 G4double yOutput[],
68 G4double yError[]) final;
69
70 inline void StepWithFinalDerivate(const G4double yInput[],
71 const G4double dydx[],
72 G4double hstep,
73 G4double yOutput[],
74 G4double yError[],
75 G4double dydxOutput[]);
76
77 inline void SetupInterpolation() {}
78
79 void Interpolate(G4double tau, G4double yOut[]) const
80 {
81 Interpolate4thOrder(yOut, tau);
82 }
83 // For calculating the output at the tau fraction of Step
84
85 G4double DistChord() const final;
86
87 inline G4int IntegratorOrder() const override { return 4; }
88
89 inline const field_utils::ShortState<N>& GetYOut() const { return fyOut; }
90
91 void Interpolate4thOrder(G4double yOut[], G4double tau) const;
92
94 void Interpolate5thOrder(G4double yOut[], G4double tau) const;
95
96 // __attribute__((always_inline))
97 inline void RightHandSideInl( const G4double y[],
98 G4double dydx[] )
99 {
100 fEquation_Rhs->T_Equation::RightHandSide(y, dydx);
101 }
102
103 inline void Stepper(const G4double yInput[],
104 const G4double dydx[],
105 G4double hstep, G4double yOutput[],
106 G4double yError[], G4double dydxOutput[])
107 {
108 StepWithFinalDerivate(yInput, dydx, hstep, yOutput, yError, dydxOutput);
109 }
110
111 T_Equation* GetSpecificEquation() { return fEquation_Rhs; }
112
113 static constexpr G4int N8 = N > 8 ? N : 8; // y[
114
115 private:
116
117 field_utils::ShortState<N> ak2, ak3, ak4, ak5, ak6, ak7, ak8, ak9;
119 field_utils::ShortState<N> fyOut, fdydxIn;
120
121 // - Simpler :
122 // field_utils::State ak2, ak3, ak4, ak5, ak6, ak7, ak8, ak9;
123 // field_utils::State fyIn, fyOut, fdydxIn;
124
125 G4double fLastStepLength = -1.0;
126 T_Equation* fEquation_Rhs;
127};
128
129// --------------------------------------------------------------------
130// G4TDormandPrince745 implementation -- borrowed from G4DormandPrince745
131//
132// DormandPrince7 - 5(4) non-FSAL
133// definition of the stepper() method that evaluates one step in
134// field propagation.
135// The coefficients and the algorithm have been adapted from
136//
137// J. R. Dormand and P. J. Prince, "A family of embedded Runge-Kutta formulae"
138// Journal of computational and applied Math., vol.6, no.1, pp.19-26, 1980.
139//
140// The Butcher table of the Dormand-Prince-7-4-5 method is as follows :
141//
142// 0 |
143// 1/5 | 1/5
144// 3/10| 3/40 9/40
145// 4/5 | 44/45 56/15 32/9
146// 8/9 | 19372/6561 25360/2187 64448/6561 212/729
147// 1 | 9017/3168 355/33 46732/5247 49/176 5103/18656
148// 1 | 35/384 0 500/1113 125/192 2187/6784 11/84
149// ------------------------------------------------------------------------
150// 35/384 0 500/1113 125/192 2187/6784 11/84 0
151// 5179/57600 0 7571/16695 393/640 92097/339200 187/2100 1/40
152//
153// --------------------------------------------------------------------
154
155// Constructor
156//
157template <class T_Equation, unsigned int N>
159 : G4MagIntegratorStepper(dynamic_cast<G4EquationOfMotion*>(equation), N )
160 , fEquation_Rhs(equation)
161{
162 // assert( dynamic_cast<G4EquationOfMotion*>(equation) != nullptr );
163 if( dynamic_cast<G4EquationOfMotion*>(equation) == nullptr )
164 {
165 G4Exception("G4TDormandPrince745: constructor", "GeomField0001",
166 FatalException, "T_Equation is not an G4EquationOfMotion.");
167 }
168
169 /***
170 assert( equation->GetNumberOfVariables == N );
171 if( equation->GetNumberOfVariables != N ){
172 G4ExceptionDescription msg;
173 msg << "Equation has an incompatible number of variables." ;
174 msg << " template N = " << N << " equation-Nvar= "
175 << equation->GetNumberOfVariables;
176 G4Exception("G4TCashKarpRKF45: constructor", "GeomField0001",
177 FatalException, msg );
178 } ****/
179}
180
181template <class T_Equation, unsigned int N>
183 G4TDormandPrince45(T_Equation* equation, G4int numVar )
184 : G4TDormandPrince45<T_Equation,N>(equation )
185{
186 if( numVar != G4int(N))
187 {
189 msg << "Equation has an incompatible number of variables." ;
190 msg << " template N = " << N
191 << " argument numVar = " << numVar ;
192 // << " equation-Nvar= " << equation->GetNumberOfVariables(); // --> Expected later
193 G4Exception("G4TCashKarpRKF45: constructor", "GeomField0001",
195 }
196 assert( numVar == N );
197}
198
199template <class T_Equation, unsigned int N>
200inline void
202 const G4double dydx[],
203 G4double hstep,
204 G4double yOutput[],
205 G4double yError[],
206 G4double dydxOutput[])
207{
208 StepWithError(yInput, dydx, hstep, yOutput, yError);
209 field_utils::copy(dydxOutput, ak7, N);
210}
211
212// Stepper
213//
214// Passing in the value of yInput[],the first time dydx[] and Step length
215// Giving back yOut and yErr arrays for output and error respectively
216//
217
218template <class T_Equation, unsigned int N>
219inline void
221 const G4double dydx[],
222 G4double hstep,
223 G4double yOut[],
224 G4double yErr[] )
225{
226 // The parameters of the Butcher tableu
227 //
228 constexpr G4double b21 = 0.2,
229 b31 = 3.0 / 40.0, b32 = 9.0 / 40.0,
230 b41 = 44.0 / 45.0, b42 = -56.0 / 15.0, b43 = 32.0/9.0,
231
232 b51 = 19372.0 / 6561.0, b52 = -25360.0 / 2187.0, b53 = 64448.0 / 6561.0,
233 b54 = -212.0 / 729.0,
234
235 b61 = 9017.0 / 3168.0 , b62 = -355.0 / 33.0,
236 b63 = 46732.0 / 5247.0, b64 = 49.0 / 176.0,
237 b65 = -5103.0 / 18656.0,
238
239 b71 = 35.0 / 384.0, b72 = 0.,
240 b73 = 500.0 / 1113.0, b74 = 125.0 / 192.0,
241 b75 = -2187.0 / 6784.0, b76 = 11.0 / 84.0,
242
243 // Sum of columns, sum(bij) = ei
244 // e1 = 0. ,
245 // e2 = 1.0/5.0 ,
246 // e3 = 3.0/10.0 ,
247 // e4 = 4.0/5.0 ,
248 // e5 = 8.0/9.0 ,
249 // e6 = 1.0 ,
250 // e7 = 1.0 ,
251
252 // Difference between the higher and the lower order method coeff. :
253 // b7j are the coefficients of higher order
254
255 dc1 = -(b71 - 5179.0 / 57600.0),
256 dc2 = -(b72 - .0),
257 dc3 = -(b73 - 7571.0 / 16695.0),
258 dc4 = -(b74 - 393.0 / 640.0),
259 dc5 = -(b75 + 92097.0 / 339200.0),
260 dc6 = -(b76 - 187.0 / 2100.0),
261 dc7 = -(- 1.0 / 40.0);
262
263 // const G4int numberOfVariables = GetNumberOfVariables();
264 // The number of variables to be integrated over
266
267 yOut[7] = yTemp[7] = fyIn[7] = yInput[7]; // Pass along the time - used in RightHandSide
268
269 // Saving yInput because yInput and yOut can be aliases for same array
270 //
271 for(unsigned int i = 0; i < N; ++i)
272 {
273 fyIn[i] = yInput[i];
274 yTemp[i] = yInput[i] + b21 * hstep * dydx[i];
275 }
276 RightHandSideInl(yTemp, ak2); // 2nd stage
277
278 for(unsigned int i = 0; i < N; ++i)
279 {
280 yTemp[i] = fyIn[i] + hstep * (b31 * dydx[i] + b32 * ak2[i]);
281 }
282 RightHandSideInl(yTemp, ak3); // 3rd stage
283
284 for(unsigned int i = 0; i < N; ++i)
285 {
286 yTemp[i] = fyIn[i] + hstep * (
287 b41 * dydx[i] + b42 * ak2[i] + b43 * ak3[i]);
288 }
289 RightHandSideInl(yTemp, ak4); // 4th stage
290
291 for(unsigned int i = 0; i < N; ++i)
292 {
293 yTemp[i] = fyIn[i] + hstep * (
294 b51 * dydx[i] + b52 * ak2[i] + b53 * ak3[i] + b54 * ak4[i]);
295 }
296 RightHandSideInl(yTemp, ak5); // 5th stage
297
298 for(unsigned int i = 0; i < N; ++i)
299 {
300 yTemp[i] = fyIn[i] + hstep * (
301 b61 * dydx[i] + b62 * ak2[i] +
302 b63 * ak3[i] + b64 * ak4[i] + b65 * ak5[i]);
303 }
304 RightHandSideInl(yTemp, ak6); // 6th stage
305
306 for(unsigned int i = 0; i < N; ++i)
307 {
308 yOut[i] = fyIn[i] + hstep * (
309 b71 * dydx[i] + b72 * ak2[i] + b73 * ak3[i] +
310 b74 * ak4[i] + b75 * ak5[i] + b76 * ak6[i]);
311 }
312 RightHandSideInl(yOut, ak7); // 7th and Final stage
313
314 for(unsigned int i = 0; i < N; ++i)
315 {
316 yErr[i] = hstep * (
317 dc1 * dydx[i] + dc2 * ak2[i] +
318 dc3 * ak3[i] + dc4 * ak4[i] +
319 dc5 * ak5[i] + dc6 * ak6[i] + dc7 * ak7[i]
320 ) + 1.5e-18;
321
322 // Store Input and Final values, for possible use in calculating chord
323 //
324 fyOut[i] = yOut[i];
325 fdydxIn[i] = dydx[i];
326 }
327
328 fLastStepLength = hstep;
329}
330
331template <class T_Equation, unsigned int N >
332inline void
334 const G4double dydx[],
335 G4double Step,
336 G4double yOutput[],
337 G4double yError[])
338{
339 assert( yOutput != yInput );
340 assert( yError != yInput );
341
342 StepWithError( yInput, dydx, Step, yOutput, yError);
343}
344
345template <class T_Equation, unsigned int N>
347{
348 // Coefficients were taken from Some Practical Runge-Kutta Formulas
349 // by Lawrence F. Shampine, page 149, c*
350 //
351 const G4double hf1 = 6025192743.0 / 30085553152.0,
352 hf3 = 51252292925.0 / 65400821598.0,
353 hf4 = - 2691868925.0 / 45128329728.0,
354 hf5 = 187940372067.0 / 1594534317056.0,
355 hf6 = - 1776094331.0 / 19743644256.0,
356 hf7 = 11237099.0 / 235043384.0;
357
358 G4ThreeVector mid;
359
360 for(unsigned int i = 0; i < 3; ++i)
361 {
362 mid[i] = fyIn[i] + 0.5 * fLastStepLength * (
363 hf1 * fdydxIn[i] + hf3 * ak3[i] +
364 hf4 * ak4[i] + hf5 * ak5[i] + hf6 * ak6[i] + hf7 * ak7[i]);
365 }
366
367 const G4ThreeVector begin = makeVector(fyIn, field_utils::Value3D::Position);
368 const G4ThreeVector end = makeVector(fyOut, field_utils::Value3D::Position);
369
370 return G4LineSection::Distline(mid, begin, end);
371}
372
373// The lower (4th) order interpolant given by Dormand and Prince:
374// J. R. Dormand and P. J. Prince, "Runge-Kutta triples"
375// Computers & Mathematics with Applications, vol. 12, no. 9,
376// pp. 1007-1017, 1986.
377//
378template <class T_Equation, unsigned int N>
379inline void
381 G4double tau) const
382{
383 const G4double tau2 = tau * tau,
384 tau3 = tau * tau2,
385 tau4 = tau2 * tau2;
386
387 const G4double bf1 = 1.0 / 11282082432.0 * (
388 157015080.0 * tau4 - 13107642775.0 * tau3 + 34969693132.0 * tau2 -
389 32272833064.0 * tau + 11282082432.0);
390
391 const G4double bf3 = - 100.0 / 32700410799.0 * tau * (
392 15701508.0 * tau3 - 914128567.0 * tau2 + 2074956840.0 * tau -
393 1323431896.0);
394
395 const G4double bf4 = 25.0 / 5641041216.0 * tau * (
396 94209048.0 * tau3 - 1518414297.0 * tau2 + 2460397220.0 * tau -
397 889289856.0);
398
399 const G4double bf5 = - 2187.0 / 199316789632.0 * tau * (
400 52338360.0 * tau3 - 451824525.0 * tau2 + 687873124.0 * tau -
401 259006536.0);
402
403 const G4double bf6 = 11.0 / 2467955532.0 * tau * (
404 106151040.0 * tau3 - 661884105.0 * tau2 +
405 946554244.0 * tau - 361440756.0);
406
407 const G4double bf7 = 1.0 / 29380423.0 * tau * (1.0 - tau) * (
408 8293050.0 * tau2 - 82437520.0 * tau + 44764047.0);
409
410 for(unsigned int i = 0; i < N; ++i)
411 {
412 yOut[i] = fyIn[i] + fLastStepLength * tau * (
413 bf1 * fdydxIn[i] + bf3 * ak3[i] + bf4 * ak4[i] +
414 bf5 * ak5[i] + bf6 * ak6[i] + bf7 * ak7[i]);
415 }
416}
417
418// Following interpolant of order 5 was given by Baker,Dormand,Gilmore, Prince :
419// T. S. Baker, J. R. Dormand, J. P. Gilmore, and P. J. Prince,
420// "Continuous approximation with embedded Runge-Kutta methods"
421// Applied Numerical Mathematics, vol. 22, no. 1, pp. 51-62, 1996.
422//
423// Calculating the extra stages for the interpolant
424//
425template <class T_Equation, unsigned int N>
427{
428 // Coefficients for the additional stages
429 //
430 const G4double b81 = 6245.0 / 62208.0,
431 b82 = 0.0,
432 b83 = 8875.0 / 103032.0,
433 b84 = -125.0 / 1728.0,
434 b85 = 801.0 / 13568.0,
435 b86 = -13519.0 / 368064.0,
436 b87 = 11105.0 / 368064.0,
437
438 b91 = 632855.0 / 4478976.0,
439 b92 = 0.0,
440 b93 = 4146875.0 / 6491016.0,
441 b94 = 5490625.0 /14183424.0,
442 b95 = -15975.0 / 108544.0,
443 b96 = 8295925.0 / 220286304.0,
444 b97 = -1779595.0 / 62938944.0,
445 b98 = -805.0 / 4104.0;
446
448
449 // Evaluate the extra stages
450 //
451 for(unsigned int i = 0; i < N; ++i)
452 {
453 yTemp[i] = fyIn[i] + fLastStepLength * (
454 b81 * fdydxIn[i] + b82 * ak2[i] + b83 * ak3[i] +
455 b84 * ak4[i] + b85 * ak5[i] + b86 * ak6[i] +
456 b87 * ak7[i]
457 );
458 }
459 RightHandSideInl(yTemp, ak8); // 8th Stage
460
461 for(unsigned int i = 0; i < N; ++i)
462 {
463 yTemp[i] = fyIn[i] + fLastStepLength * (
464 b91 * fdydxIn[i] + b92 * ak2[i] + b93 * ak3[i] +
465 b94 * ak4[i] + b95 * ak5[i] + b96 * ak6[i] +
466 b97 * ak7[i] + b98 * ak8[i]
467 );
468 }
469 RightHandSideInl(yTemp, ak9); // 9th Stage
470}
471
472// Calculating the interpolated result yOut with the coefficients
473//
474template <class T_Equation, unsigned int N>
476Interpolate5thOrder(G4double yOut[], G4double tau) const
477{
478 // Define the coefficients for the polynomials
479 //
480 G4double bi[10][5];
481
482 // COEFFICIENTS OF bi[1]
483 bi[1][0] = 1.0,
484 bi[1][1] = -38039.0 / 7040.0,
485 bi[1][2] = 125923.0 / 10560.0,
486 bi[1][3] = -19683.0 / 1760.0,
487 bi[1][4] = 3303.0 / 880.0,
488 // --------------------------------------------------------
489 //
490 // COEFFICIENTS OF bi[2]
491 bi[2][0] = 0.0,
492 bi[2][1] = 0.0,
493 bi[2][2] = 0.0,
494 bi[2][3] = 0.0,
495 bi[2][4] = 0.0,
496 // --------------------------------------------------------
497 //
498 // COEFFICIENTS OF bi[3]
499 bi[3][0] = 0.0,
500 bi[3][1] = -12500.0 / 4081.0,
501 bi[3][2] = 205000.0 / 12243.0,
502 bi[3][3] = -90000.0 / 4081.0,
503 bi[3][4] = 36000.0 / 4081.0,
504 // --------------------------------------------------------
505 //
506 // COEFFICIENTS OF bi[4]
507 bi[4][0] = 0.0,
508 bi[4][1] = -3125.0 / 704.0,
509 bi[4][2] = 25625.0 / 1056.0,
510 bi[4][3] = -5625.0 / 176.0,
511 bi[4][4] = 1125.0 / 88.0,
512 // --------------------------------------------------------
513 //
514 // COEFFICIENTS OF bi[5]
515 bi[5][0] = 0.0,
516 bi[5][1] = 164025.0 / 74624.0,
517 bi[5][2] = -448335.0 / 37312.0,
518 bi[5][3] = 295245.0 / 18656.0,
519 bi[5][4] = -59049.0 / 9328.0,
520 // --------------------------------------------------------
521 //
522 // COEFFICIENTS OF bi[6]
523 bi[6][0] = 0.0,
524 bi[6][1] = -25.0 / 28.0,
525 bi[6][2] = 205.0 / 42.0,
526 bi[6][3] = -45.0 / 7.0,
527 bi[6][4] = 18.0 / 7.0,
528 // --------------------------------------------------------
529 //
530 // COEFFICIENTS OF bi[7]
531 bi[7][0] = 0.0,
532 bi[7][1] = -2.0 / 11.0,
533 bi[7][2] = 73.0 / 55.0,
534 bi[7][3] = -171.0 / 55.0,
535 bi[7][4] = 108.0 / 55.0,
536 // --------------------------------------------------------
537 //
538 // COEFFICIENTS OF bi[8]
539 bi[8][0] = 0.0,
540 bi[8][1] = 189.0 / 22.0,
541 bi[8][2] = -1593.0 / 55.0,
542 bi[8][3] = 3537.0 / 110.0,
543 bi[8][4] = -648.0 / 55.0,
544 // --------------------------------------------------------
545 //
546 // COEFFICIENTS OF bi[9]
547 bi[9][0] = 0.0,
548 bi[9][1] = 351.0 / 110.0,
549 bi[9][2] = -999.0 / 55.0,
550 bi[9][3] = 2943.0 / 110.0,
551 bi[9][4] = -648.0 / 55.0;
552 // --------------------------------------------------------
553
554 // Calculating the polynomials
555
556 G4double b[10];
557 std::memset(b, 0.0, sizeof(b));
558
559 G4double tauPower = 1.0;
560 for(G4int j = 0; j <= 4; ++j)
561 {
562 for(G4int iStage = 1; iStage <= 9; ++iStage)
563 {
564 b[iStage] += bi[iStage][j] * tauPower;
565 }
566 tauPower *= tau;
567 }
568
569 const G4double stepLen = fLastStepLength * tau;
570 for(G4int i = 0; i < N; ++i)
571 {
572 yOut[i] = fyIn[i] + stepLen * (
573 b[1] * fdydxIn[i] + b[2] * ak2[i] + b[3] * ak3[i] +
574 b[4] * ak4[i] + b[5] * ak5[i] + b[6] * ak6[i] +
575 b[7] * ak7[i] + b[8] * ak8[i] + b[9] * ak9[i]
576 );
577 }
578}
579
580#endif
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
G4int IntegratorOrder() const override
void RightHandSideInl(const G4double y[], G4double dydx[])
G4TDormandPrince45(T_Equation *equation)
T_Equation * GetSpecificEquation()
void Interpolate(G4double tau, G4double yOut[]) const
G4double DistChord() const final
void Interpolate4thOrder(G4double yOut[], G4double tau) const
void Stepper(const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[]) final
void StepWithFinalDerivate(const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[], G4double dydxOutput[])
void Stepper(const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[], G4double dydxOutput[])
const field_utils::ShortState< N > & GetYOut() const
static constexpr G4int N8
void Interpolate5thOrder(G4double yOut[], G4double tau) const
void StepWithError(const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[])
#define N
Definition crc32.c:57
#define inline
Definition internal.h:104
void copy(G4double dst[], const G4double src[], std::size_t size=G4FieldTrack::ncompSVEC)
G4double[N] ShortState