Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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
43#include <cassert>
45#include "G4FieldUtils.hh"
46
47template <class T_Equation, unsigned int N = 6 >
49{
50 public:
51
52 G4TDormandPrince45(T_Equation* equation );
53 G4TDormandPrince45(T_Equation* equation, G4int numVar ); // must have numVar == N
54
55 inline void
56 StepWithError(const G4double yInput[],
57 const G4double dydx[],
58 G4double hstep,
59 G4double yOutput[],
60 G4double yError[] ) ;
61
62 virtual void Stepper(const G4double yInput[],
63 const G4double dydx[],
64 G4double hstep,
65 G4double yOutput[],
66 G4double yError[]) override final;
67
68 inline
69 void StepWithFinalDerivate(const G4double yInput[],
70 const G4double dydx[],
71 G4double hstep,
72 G4double yOutput[],
73 G4double yError[],
74 G4double dydxOutput[]);
75
76 inline void SetupInterpolation() {}
77
78 void Interpolate(G4double tau, G4double yOut[]) const
79 {
80 Interpolate4thOrder(yOut, tau);
81 }
82 // For calculating the output at the tau fraction of Step
83
84 virtual G4double DistChord() const override final;
85
86 virtual G4int IntegratorOrder() const override { return 4; }
87
88 const field_utils::ShortState<N>& GetYOut() const { return fyOut; }
89
90 void Interpolate4thOrder(G4double yOut[], G4double tau) const;
91
93 void Interpolate5thOrder(G4double yOut[], G4double tau) const;
94
95 // __attribute__((always_inline))
96 void RightHandSideInl( const G4double y[],
97 G4double dydx[] )
98 {
99 fEquation_Rhs->T_Equation::RightHandSide(y, dydx);
100 }
101
102 inline
103 void Stepper(const G4double yInput[], const G4double dydx[],
104 G4double hstep, G4double yOutput[],
105 G4double yError[], G4double dydxOutput[])
106 {
107 StepWithFinalDerivate(yInput, dydx, hstep,
108 yOutput, yError, dydxOutput);
109 }
110
111 T_Equation* GetSpecificEquation() { return fEquation_Rhs; }
112
113 static constexpr int 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// G4TDormandPrince745 implementation -- borrowed from G4DormandPrince745
130//
131// DormandPrince7 - 5(4) non-FSAL
132// definition of the stepper() method that evaluates one step in
133// field propagation.
134// The coefficients and the algorithm have been adapted from
135//
136// J. R. Dormand and P. J. Prince, "A family of embedded Runge-Kutta formulae"
137// Journal of computational and applied Math., vol.6, no.1, pp.19-26, 1980.
138//
139// The Butcher table of the Dormand-Prince-7-4-5 method is as follows :
140//
141// 0 |
142// 1/5 | 1/5
143// 3/10| 3/40 9/40
144// 4/5 | 44/45 56/15 32/9
145// 8/9 | 19372/6561 25360/2187 64448/6561 212/729
146// 1 | 9017/3168 355/33 46732/5247 49/176 5103/18656
147// 1 | 35/384 0 500/1113 125/192 2187/6784 11/84
148// ------------------------------------------------------------------------
149// 35/384 0 500/1113 125/192 2187/6784 11/84 0
150// 5179/57600 0 7571/16695 393/640 92097/339200 187/2100 1/40
151//
152// Created: Somnath Banerjee, Google Summer of Code 2015, 25 May 2015
153// Supervision: John Apostolakis, CERN
154// --------------------------------------------------------------------
155
156#include "G4LineSection.hh"
157
158#include <cstring>
159
160// using namespace field_utils;
161
162/////////////////////////////////////////////////////////////////////
163// Constructor
164//
165template <class T_Equation, unsigned int N>
167 : G4MagIntegratorStepper(dynamic_cast<G4EquationOfMotion*>(equation), N )
168 , fEquation_Rhs(equation)
169{
170 // assert( dynamic_cast<G4EquationOfMotion*>(equation) != nullptr );
171 if( dynamic_cast<G4EquationOfMotion*>(equation) == nullptr )
172 {
173 G4Exception("G4TDormandPrince745: constructor", "GeomField0001",
174 FatalException, "T_Equation is not an G4EquationOfMotion.");
175 }
176
177 /***
178 assert( equation->GetNumberOfVariables == N );
179 if( equation->GetNumberOfVariables != N ){
180 G4ExceptionDescription msg;
181 msg << "Equation has an incompatible number of variables." ;
182 msg << " template N = " << N << " equation-Nvar= "
183 << equation->GetNumberOfVariables;
184 G4Exception("G4TCashKarpRKF45: constructor", "GeomField0001",
185 FatalException, msg );
186 } ****/
187}
188
189template <class T_Equation, unsigned int N>
191 : G4TDormandPrince45<T_Equation,N>(equation )
192{
193 if( numVar != G4int(N)){
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>
207StepWithFinalDerivate(const G4double yInput[],
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
337
338template <class T_Equation, unsigned int N >
339inline void
341 const G4double dydx[],
342 G4double Step,
343 G4double yOutput[],
344 G4double yError[])
345{
346 assert( yOutput != yInput );
347 assert( yError != yInput );
348
349 StepWithError( yInput, dydx, Step, yOutput, yError);
350}
351
352template <class T_Equation, unsigned int N>
355{
356 // Coefficients were taken from Some Practical Runge-Kutta Formulas
357 // by Lawrence F. Shampine, page 149, c*
358 //
359 const G4double hf1 = 6025192743.0 / 30085553152.0,
360 hf3 = 51252292925.0 / 65400821598.0,
361 hf4 = - 2691868925.0 / 45128329728.0,
362 hf5 = 187940372067.0 / 1594534317056.0,
363 hf6 = - 1776094331.0 / 19743644256.0,
364 hf7 = 11237099.0 / 235043384.0;
365
366 G4ThreeVector mid;
367
368 for(unsigned int i = 0; i < 3; ++i)
369 {
370 mid[i] = fyIn[i] + 0.5 * fLastStepLength * (
371 hf1 * fdydxIn[i] + hf3 * ak3[i] +
372 hf4 * ak4[i] + hf5 * ak5[i] + hf6 * ak6[i] + hf7 * ak7[i]);
373 }
374
375 const G4ThreeVector begin = makeVector(fyIn, field_utils::Value3D::Position);
376 const G4ThreeVector end = makeVector(fyOut, field_utils::Value3D::Position);
377
378 return G4LineSection::Distline(mid, begin, end);
379}
380
381// The lower (4th) order interpolant given by Dormand and Prince:
382// J. R. Dormand and P. J. Prince, "Runge-Kutta triples"
383// Computers & Mathematics with Applications, vol. 12, no. 9,
384// pp. 1007-1017, 1986.
385//
386template <class T_Equation, unsigned int N>
387void
389Interpolate4thOrder(G4double yOut[], G4double tau) const
390{
391 // const G4int numberOfVariables = GetNumberOfVariables();
392
393 const G4double tau2 = tau * tau,
394 tau3 = tau * tau2,
395 tau4 = tau2 * tau2;
396
397 const G4double bf1 = 1.0 / 11282082432.0 * (
398 157015080.0 * tau4 - 13107642775.0 * tau3 + 34969693132.0 * tau2 -
399 32272833064.0 * tau + 11282082432.0);
400
401 const G4double bf3 = - 100.0 / 32700410799.0 * tau * (
402 15701508.0 * tau3 - 914128567.0 * tau2 + 2074956840.0 * tau -
403 1323431896.0);
404
405 const G4double bf4 = 25.0 / 5641041216.0 * tau * (
406 94209048.0 * tau3 - 1518414297.0 * tau2 + 2460397220.0 * tau -
407 889289856.0);
408
409 const G4double bf5 = - 2187.0 / 199316789632.0 * tau * (
410 52338360.0 * tau3 - 451824525.0 * tau2 + 687873124.0 * tau -
411 259006536.0);
412
413 const G4double bf6 = 11.0 / 2467955532.0 * tau * (
414 106151040.0 * tau3 - 661884105.0 * tau2 +
415 946554244.0 * tau - 361440756.0);
416
417 const G4double bf7 = 1.0 / 29380423.0 * tau * (1.0 - tau) * (
418 8293050.0 * tau2 - 82437520.0 * tau + 44764047.0);
419
420 for(unsigned int i = 0; i < N; ++i)
421 {
422 yOut[i] = fyIn[i] + fLastStepLength * tau * (
423 bf1 * fdydxIn[i] + bf3 * ak3[i] + bf4 * ak4[i] +
424 bf5 * ak5[i] + bf6 * ak6[i] + bf7 * ak7[i]);
425 }
426}
427
428// Following interpolant of order 5 was given by Baker,Dormand,Gilmore, Prince :
429// T. S. Baker, J. R. Dormand, J. P. Gilmore, and P. J. Prince,
430// "Continuous approximation with embedded Runge-Kutta methods"
431// Applied Numerical Mathematics, vol. 22, no. 1, pp. 51-62, 1996.
432//
433// Calculating the extra stages for the interpolant
434//
435template <class T_Equation, unsigned int N>
437{
438 // Coefficients for the additional stages
439 //
440 const G4double b81 = 6245.0 / 62208.0,
441 b82 = 0.0,
442 b83 = 8875.0 / 103032.0,
443 b84 = -125.0 / 1728.0,
444 b85 = 801.0 / 13568.0,
445 b86 = -13519.0 / 368064.0,
446 b87 = 11105.0 / 368064.0,
447
448 b91 = 632855.0 / 4478976.0,
449 b92 = 0.0,
450 b93 = 4146875.0 / 6491016.0,
451 b94 = 5490625.0 /14183424.0,
452 b95 = -15975.0 / 108544.0,
453 b96 = 8295925.0 / 220286304.0,
454 b97 = -1779595.0 / 62938944.0,
455 b98 = -805.0 / 4104.0;
456
457 // const G4int numberOfVariables = GetNumberOfVariables();
459
460 // Evaluate the extra stages
461 //
462 for(unsigned int i = 0; i < N; ++i)
463 {
464 yTemp[i] = fyIn[i] + fLastStepLength * (
465 b81 * fdydxIn[i] + b82 * ak2[i] + b83 * ak3[i] +
466 b84 * ak4[i] + b85 * ak5[i] + b86 * ak6[i] +
467 b87 * ak7[i]
468 );
469 }
470 RightHandSideInl(yTemp, ak8); // 8th Stage
471
472 for(unsigned int i = 0; i < N; ++i)
473 {
474 yTemp[i] = fyIn[i] + fLastStepLength * (
475 b91 * fdydxIn[i] + b92 * ak2[i] + b93 * ak3[i] +
476 b94 * ak4[i] + b95 * ak5[i] + b96 * ak6[i] +
477 b97 * ak7[i] + b98 * ak8[i]
478 );
479 }
480 RightHandSideInl(yTemp, ak9); // 9th Stage
481}
482
483// Calculating the interpolated result yOut with the coefficients
484//
485template <class T_Equation, unsigned int N>
487Interpolate5thOrder(G4double yOut[], G4double tau) const
488{
489 // Define the coefficients for the polynomials
490 //
491 G4double bi[10][5];
492
493 // COEFFICIENTS OF bi[1]
494 bi[1][0] = 1.0,
495 bi[1][1] = -38039.0 / 7040.0,
496 bi[1][2] = 125923.0 / 10560.0,
497 bi[1][3] = -19683.0 / 1760.0,
498 bi[1][4] = 3303.0 / 880.0,
499 // --------------------------------------------------------
500 //
501 // COEFFICIENTS OF bi[2]
502 bi[2][0] = 0.0,
503 bi[2][1] = 0.0,
504 bi[2][2] = 0.0,
505 bi[2][3] = 0.0,
506 bi[2][4] = 0.0,
507 // --------------------------------------------------------
508 //
509 // COEFFICIENTS OF bi[3]
510 bi[3][0] = 0.0,
511 bi[3][1] = -12500.0 / 4081.0,
512 bi[3][2] = 205000.0 / 12243.0,
513 bi[3][3] = -90000.0 / 4081.0,
514 bi[3][4] = 36000.0 / 4081.0,
515 // --------------------------------------------------------
516 //
517 // COEFFICIENTS OF bi[4]
518 bi[4][0] = 0.0,
519 bi[4][1] = -3125.0 / 704.0,
520 bi[4][2] = 25625.0 / 1056.0,
521 bi[4][3] = -5625.0 / 176.0,
522 bi[4][4] = 1125.0 / 88.0,
523 // --------------------------------------------------------
524 //
525 // COEFFICIENTS OF bi[5]
526 bi[5][0] = 0.0,
527 bi[5][1] = 164025.0 / 74624.0,
528 bi[5][2] = -448335.0 / 37312.0,
529 bi[5][3] = 295245.0 / 18656.0,
530 bi[5][4] = -59049.0 / 9328.0,
531 // --------------------------------------------------------
532 //
533 // COEFFICIENTS OF bi[6]
534 bi[6][0] = 0.0,
535 bi[6][1] = -25.0 / 28.0,
536 bi[6][2] = 205.0 / 42.0,
537 bi[6][3] = -45.0 / 7.0,
538 bi[6][4] = 18.0 / 7.0,
539 // --------------------------------------------------------
540 //
541 // COEFFICIENTS OF bi[7]
542 bi[7][0] = 0.0,
543 bi[7][1] = -2.0 / 11.0,
544 bi[7][2] = 73.0 / 55.0,
545 bi[7][3] = -171.0 / 55.0,
546 bi[7][4] = 108.0 / 55.0,
547 // --------------------------------------------------------
548 //
549 // COEFFICIENTS OF bi[8]
550 bi[8][0] = 0.0,
551 bi[8][1] = 189.0 / 22.0,
552 bi[8][2] = -1593.0 / 55.0,
553 bi[8][3] = 3537.0 / 110.0,
554 bi[8][4] = -648.0 / 55.0,
555 // --------------------------------------------------------
556 //
557 // COEFFICIENTS OF bi[9]
558 bi[9][0] = 0.0,
559 bi[9][1] = 351.0 / 110.0,
560 bi[9][2] = -999.0 / 55.0,
561 bi[9][3] = 2943.0 / 110.0,
562 bi[9][4] = -648.0 / 55.0;
563 // --------------------------------------------------------
564
565 // Calculating the polynomials
566
567 G4double b[10];
568 std::memset(b, 0.0, sizeof(b));
569
570 G4double tauPower = 1.0;
571 for(G4int j = 0; j <= 4; ++j)
572 {
573 for(G4int iStage = 1; iStage <= 9; ++iStage)
574 {
575 b[iStage] += bi[iStage][j] * tauPower;
576 }
577 tauPower *= tau;
578 }
579
580 // const G4int numberOfVariables = GetNumberOfVariables();
581 const G4double stepLen = fLastStepLength * tau;
582 for(G4int i = 0; i < N; ++i)
583 {
584 yOut[i] = fyIn[i] + stepLen * (
585 b[1] * fdydxIn[i] + b[2] * ak2[i] + b[3] * ak3[i] +
586 b[4] * ak4[i] + b[5] * ak5[i] + b[6] * ak6[i] +
587 b[7] * ak7[i] + b[8] * ak8[i] + b[9] * ak9[i]
588 );
589 }
590}
591
592
593#endif
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
void RightHandSideInl(const G4double y[], G4double dydx[])
G4TDormandPrince45(T_Equation *equation)
T_Equation * GetSpecificEquation()
void Interpolate(G4double tau, G4double yOut[]) const
virtual G4int IntegratorOrder() const override
static constexpr int N8
void Interpolate4thOrder(G4double yOut[], G4double tau) const
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
virtual G4double DistChord() const override final
virtual void Stepper(const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[]) override final
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:56
void copy(G4double dst[], const G4double src[], std::size_t size=G4FieldTrack::ncompSVEC)
Definition: G4FieldUtils.cc:98
G4double[N] ShortState
Definition: G4FieldUtils.hh:48