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