Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4FSALDormandPrince745.cc
이 파일의 문서화 페이지로 가기
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 // DormandPrince7 - 5(4) implementation by Somnath Banerjee
27 // Supervision / code review: John Apostolakis
28 //
29 // Sponsored by Google in Google Summer of Code 2015.
30 //
31 // First version: 25 May 2015
32 //
33 // G4FSALDormandPrince745.cc
34 // Geant4
35 //
36 // This is the source file of G4FSALDormandPrince745 class containing the
37 // definition of the stepper() method that evaluates one step in
38 // field propagation.
39 // The Butcher table of the FDormand-Prince-7-4-5 method is as follows :
40 //
41 // 0 |
42 // 1/5 | 1/5
43 // 3/10| 3/40 9/40
44 // 4/5 | 44/45 −56/15 32/9
45 // 8/9 | 19372/6561 −25360/2187 64448/6561 −212/729
46 // 1 | 9017/3168 −355/33 46732/5247 49/176 −5103/18656
47 // 1 | 35/384 0 500/1113 125/192 −2187/6784 11/84
48 // ---------------------------------------------------------------------------
49 // 35/384 0 500/1113 125/192 −2187/6784 11/84 0
50 // 5179/57600 0 7571/16695 393/640 −92097/339200 187/2100 1/40
51 //
52 // Implementation by Somnath Banerjee - GSoC 2015
53 // Work supported by Google as part of Google Summer of Code 2015.
54 // Supervision / code review: John Apostolakis
55 //
56 // First version: June 2015 - Somnath Banerjee
57 
59 #include "G4LineSection.hh"
60 #include <cmath>
61 
62 //Constructor
64  G4int noIntegrationVariables,
65  G4bool primary)
66  : G4VFSALIntegrationStepper(EqRhs, noIntegrationVariables)
67 {
68 
69  const G4int numberOfVariables = noIntegrationVariables;
70 
71  //New Chunk of memory being created for use by the stepper
72 
73  //aki - for storing intermediate RHS
74  ak2 = new G4double[numberOfVariables];
75  ak3 = new G4double[numberOfVariables];
76  ak4 = new G4double[numberOfVariables];
77  ak5 = new G4double[numberOfVariables];
78  ak6 = new G4double[numberOfVariables];
79  ak7 = new G4double[numberOfVariables];
80  // Also always allocate arrays for interpolation stages
81  ak8 = new G4double[numberOfVariables];
82  ak9 = new G4double[numberOfVariables];
83 
84  yTemp = new G4double[numberOfVariables] ;
85  yIn = new G4double[numberOfVariables] ;
86 
87  pseudoDydx_for_DistChord = new G4double[numberOfVariables];
88 
89  fInitialDyDx = new G4double[numberOfVariables];
90  fLastInitialVector = new G4double[numberOfVariables] ;
91  fLastFinalVector = new G4double[numberOfVariables] ;
92  fLastDyDx = new G4double[numberOfVariables];
93 
94  fMidVector = new G4double[numberOfVariables];
95  fMidError = new G4double[numberOfVariables];
96 
97  fAuxStepper = nullptr;
98  if( primary )
99  {
100  fAuxStepper = new G4FSALDormandPrince745(EqRhs, numberOfVariables,
101  !primary);
102  }
103  fLastStepLength = -1.0;
104 }
105 
106 
107 //Destructor
109 {
110  //clear all previously allocated memory for stepper and DistChord
111  delete[] ak2; ak2=nullptr;
112  delete[] ak3; ak3=nullptr;
113  delete[] ak4; ak4=nullptr;
114  delete[] ak5; ak5=nullptr;
115  delete[] ak6; ak6=nullptr;
116  delete[] ak7; ak7=nullptr;
117  delete[] ak8; ak8=nullptr;
118  delete[] ak9; ak9=nullptr;
119 
120  delete[] yTemp; yTemp= nullptr;
121  delete[] yIn; yIn= nullptr;
122 
124  delete[] fInitialDyDx; fInitialDyDx= nullptr;
125 
126  delete[] fLastInitialVector; fLastInitialVector= nullptr;
127  delete[] fLastFinalVector; fLastFinalVector = nullptr;
128  delete[] fLastDyDx; fLastDyDx = nullptr;
129  delete[] fMidVector; fMidVector = nullptr;
130  delete[] fMidError; fMidError = nullptr;
131 
132  delete fAuxStepper; fAuxStepper= nullptr;
133 }
134 
135 
136 //Stepper :
137 
138 // Passing in the value of yInput[],the first time dydx[] and Step length
139 // Giving back yOut and yErr arrays for output and error respectively
140 
142  const G4double dydx[],
143  G4double Step,
144  G4double yOut[],
145  G4double yErr[],
146  G4double nextDydx[]
147  )
148 {
149  G4int i;
150 
151  //The various constants defined on the basis of butcher tableu
152  const G4double //G4double - only once
153  b21 = 0.2 ,
154 
155  b31 = 3.0/40.0, b32 = 9.0/40.0 ,
156 
157  b41 = 44.0/45.0, b42 = -56.0/15.0, b43 = 32.0/9.0,
158 
159  b51 = 19372.0/6561.0, b52 = -25360.0/2187.0, b53 = 64448.0/6561.0,
160  b54 = -212.0/729.0 ,
161 
162  b61 = 9017.0/3168.0 , b62 = -355.0/33.0,
163  b63 = 46732.0/5247.0 , b64 = 49.0/176.0 ,
164  b65 = -5103.0/18656.0 ,
165 
166  b71 = 35.0/384.0, b72 = 0.,
167  b73 = 500.0/1113.0, b74 = 125.0/192.0,
168  b75 = -2187.0/6784.0, b76 = 11.0/84.0,
169 
170 // c1 = 35.0/384.0, c2 = .0,
171 // c3 = 500.0/1113.0, c4 = 125.0/192.0,
172 // c5 = -2187.0/6784.0, c6 = 11.0/84.0,
173 // c7 = 0,
174 
175  dc1 = b71 - 5179.0/57600.0,
176  dc2 = b72 - .0,
177  dc3 = b73 - 7571.0/16695.0,
178  dc4 = b74 - 393.0/640.0,
179  dc5 = b75 + 92097.0/339200.0,
180  dc6 = b76 - 187.0/2100.0,
181  dc7 = - 1.0/40.0 ; //end of declaration
182 
183 
184  const G4int numberOfVariables= this->GetNumberOfVariables();
185  // The number of variables to be integrated over
186 
187  // Saving yInput because yInput and yOut can be aliases for same array
188  for(i=0;i<numberOfVariables;i++)
189  {
190  yIn[i] = yInput[i];
191  fInitialDyDx[i] = dydx[i];
192  }
193  // Ensure that time is initialised - in case it is not integrated
194  yOut[7] = yTemp[7] = yInput[7];
195 
196  // RightHandSide(yIn, DyDx) ;
197  // 1st Step - Not doing, getting passed
198 
199  for(i=0;i<numberOfVariables;i++)
200  {
201  yTemp[i] = yIn[i] + b21*Step*fInitialDyDx[i] ;
202  }
203  RightHandSide(yTemp, ak2) ; // 2nd Step
204 
205  for(i=0;i<numberOfVariables;i++)
206  {
207  yTemp[i] = yIn[i] + Step*(b31*fInitialDyDx[i] + b32*ak2[i]) ;
208  }
209  RightHandSide(yTemp, ak3) ; // 3rd Step
210 
211  for(i=0;i<numberOfVariables;i++)
212  {
213  yTemp[i] = yIn[i] + Step*(b41*fInitialDyDx[i] + b42*ak2[i] + b43*ak3[i]) ;
214  }
215  RightHandSide(yTemp, ak4) ; // 4th Step
216 
217  for(i=0;i<numberOfVariables;i++)
218  {
219  yTemp[i] = yIn[i] + Step*(b51*fInitialDyDx[i] + b52*ak2[i] + b53*ak3[i] +
220  b54*ak4[i]) ;
221  }
222  RightHandSide(yTemp, ak5) ; // 5th Step
223 
224  for(i=0;i<numberOfVariables;i++)
225  {
226  yTemp[i] = yIn[i] + Step*(b61*fInitialDyDx[i] + b62*ak2[i] + b63*ak3[i] +
227  b64*ak4[i] + b65*ak5[i]) ;
228  }
229  RightHandSide(yTemp, ak6) ; // 6th Step
230 
231  for(i=0;i<numberOfVariables;i++)
232  {
233  yOut[i] = yIn[i] + Step*(b71*fInitialDyDx[i] + b72*ak2[i] + b73*ak3[i] +
234  b74*ak4[i] + b75*ak5[i] + b76*ak6[i]);
235  }
236  RightHandSide(yOut, ak7); //7th and Final step
237 
238  for(i=0;i<numberOfVariables;i++)
239  {
240 
241  yErr[i] = Step*(dc1*fInitialDyDx[i] + dc2*ak2[i] + dc3*ak3[i] + dc4*ak4[i] +
242  dc5*ak5[i] + dc6*ak6[i] + dc7*ak7[i] ) ;
243 
244 
245  // Store Input and Final values, for possible use in calculating chord
246  fLastInitialVector[i] = yIn[i] ;
247  fLastFinalVector[i] = yOut[i];
248  fLastDyDx[i] = fInitialDyDx[i];
249  nextDydx[i] = ak7[i];
250 
251 
252  }
253 
254  fLastStepLength = Step;
255 
256  return ;
257 }
258 
259 
260 //The following has not been tested
261 
262 //The DistChord() function fot the class - must define it here.
264 {
265  G4double distLine, distChord;
266  G4ThreeVector initialPoint, finalPoint, midPoint;
267 
268  // Store last initial and final points (they will be overwritten in self-Stepper call!)
269  initialPoint = G4ThreeVector( fLastInitialVector[0],
271  finalPoint = G4ThreeVector( fLastFinalVector[0],
273 
274  // Do half a step using StepNoErr
275 
278 
279  midPoint = G4ThreeVector( fMidVector[0], fMidVector[1], fMidVector[2]);
280 
281  // Use stored values of Initial and Endpoint + new Midpoint to evaluate
282  // distance of Chord
283 
284 
285  if (initialPoint != finalPoint)
286  {
287  distLine = G4LineSection::Distline( midPoint, initialPoint, finalPoint );
288  distChord = distLine;
289  }
290  else
291  {
292  distChord = (midPoint-initialPoint).mag();
293  }
294  return distChord;
295 }
296 
297 
299  const G4double dydx[],
300  G4double yOut[],
301  G4double Step,
302  G4double tau){
303 
304  G4double
305  bf1, bf2, bf3, bf4, bf5, bf6, bf7;
306 
307 
308 
309 
310 
311  const G4int numberOfVariables= this->GetNumberOfVariables();
312 
313  G4double tau0 = tau;
314 
315  for(int i=0;i<numberOfVariables;i++)
316  {
317  yIn[i]=yInput[i];
318  }
319 
320  G4double
321  tau_2 = tau0*tau0 ,
322  tau_3 = tau0*tau_2,
323  tau_4 = tau_2*tau_2;
324 
325  bf1 = (157015080.0*tau_4 - 13107642775.0*tau_3+ 34969693132.0*tau_2- 32272833064.0*tau
326  + 11282082432.0)/11282082432.0,
327  bf2 = 0.0 ,
328  bf3 = - 100.0*tau*(15701508.0*tau_3 - 914128567.0*tau_2 + 2074956840.0*tau
329  - 1323431896.0)/32700410799.0,
330  bf4 = 25.0*tau*(94209048.0*tau_3- 1518414297.0*tau_2+ 2460397220.0*tau - 889289856.0)/5641041216.0 ,
331  bf5 = -2187.0*tau*(52338360.0*tau_3 - 451824525.0*tau_2 + 687873124.0*tau - 259006536.0)/199316789632.0 ,
332  bf6 = 11.0*tau*(106151040.0*tau_3- 661884105.0*tau_2 + 946554244.0*tau - 361440756.0)/2467955532.0 ,
333  bf7 = tau*(1.0 - tau)*(8293050.0*tau_2 - 82437520.0*tau + 44764047.0)/ 29380423.0 ;
334 
335 
336  for( int i=0; i<numberOfVariables; i++){
337  yOut[i] = yIn[i] + Step*tau*(bf1*dydx[i] + bf2*ak2[i] + bf3*ak3[i] + bf4*ak4[i]
338  + bf5*ak5[i] + bf6*ak6[i] + bf7*ak7[i] ) ;
339  }
340 
341 
342 
343 }
344 
346  const G4double dydx[],
347  const G4double Step ){
348 
349  //Coefficients for the additional stages :
350  G4double
351  b81 = 6245.0/62208.0 ,
352  b82 = 0.0 ,
353  b83 = 8875.0/103032.0 ,
354  b84 = -125.0/1728.0 ,
355  b85 = 801.0/13568.0 ,
356  b86 = -13519.0/368064.0 ,
357  b87 = 11105.0/368064.0 ,
358 
359  b91 = 632855.0/4478976.0 ,
360  b92 = 0.0 ,
361  b93 = 4146875.0/6491016.0 ,
362  b94 = 5490625.0/14183424.0 ,
363  b95 = -15975.0/108544.0 ,
364  b96 = 8295925.0/220286304.0 ,
365  b97 = -1779595.0/62938944.0 ,
366  b98 = -805.0/4104.0 ;
367 
368  const G4int numberOfVariables= this->GetNumberOfVariables();
369 
370  // Saving yInput because yInput and yOut can be aliases for same array
371  for(int i=0;i<numberOfVariables;i++)
372  {
373  yIn[i]=yInput[i];
374  }
375 
376  yTemp[7] = yIn[7];
377 
378  //Evaluate the extra stages :
379  for(int i=0;i<numberOfVariables;i++)
380  {
381  yTemp[i] = yIn[i] + Step*( b81*dydx[i] + b82*ak2[i] + b83*ak3[i] +
382  b84*ak4[i] + b85*ak5[i] + b86*ak6[i] +
383  b87*ak7[i] );
384  }
385  RightHandSide( yTemp, ak8 ); //8th Stage
386 
387  for(int i=0;i<numberOfVariables;i++)
388  {
389  yTemp[i] = yIn[i] + Step * ( b91*dydx[i] + b92*ak2[i] + b93*ak3[i] +
390  b94*ak4[i] + b95*ak5[i] + b96*ak6[i] +
391  b97*ak7[i] + b98*ak8[i] );
392  }
393  RightHandSide( yTemp, ak9 ); //9th Stage
394 
395 
396 
397 }
398 
399 
400 
402  const G4double dydx[],
403  const G4double Step,
404  G4double yOut[],
405  G4double tau ){
406  //Define the coefficients for the polynomials
407  G4double bi[10][5], b[10];
408  G4int numberOfVariables = this->GetNumberOfVariables();
409 
410  // COEFFICIENTS OF bi[1]
411  bi[1][0] = 1.0 ,
412  bi[1][1] = -38039.0/7040.0 ,
413  bi[1][2] = 125923.0/10560.0 ,
414  bi[1][3] = -19683.0/1760.0 ,
415  bi[1][4] = 3303.0/880.0 ,
416  // --------------------------------------------------------
417  //
418  // COEFFICIENTS OF bi[2]
419  bi[2][0] = 0.0 ,
420  bi[2][1] = 0.0 ,
421  bi[2][2] = 0.0 ,
422  bi[2][3] = 0.0 ,
423  bi[2][4] = 0.0 ,
424  // --------------------------------------------------------
425  //
426  // COEFFICIENTS OF bi[3]
427  bi[3][0] = 0.0 ,
428  bi[3][1] = -12500.0/4081.0 ,
429  bi[3][2] = 205000.0/12243.0 ,
430  bi[3][3] = -90000.0/4081.0 ,
431  bi[3][4] = 36000.0/4081.0 ,
432  // --------------------------------------------------------
433  //
434  // COEFFICIENTS OF bi[4]
435  bi[4][0] = 0.0 ,
436  bi[4][1] = -3125.0/704.0 ,
437  bi[4][2] = 25625.0/1056.0 ,
438  bi[4][3] = -5625.0/176.0 ,
439  bi[4][4] = 1125.0/88.0 ,
440  // --------------------------------------------------------
441  //
442  // COEFFICIENTS OF bi[5]
443  bi[5][0] = 0.0 ,
444  bi[5][1] = 164025.0/74624.0 ,
445  bi[5][2] = -448335.0/37312.0 ,
446  bi[5][3] = 295245.0/18656.0 ,
447  bi[5][4] = -59049.0/9328.0 ,
448  // --------------------------------------------------------
449  //
450  // COEFFICIENTS OF bi[6]
451  bi[6][0] = 0.0 ,
452  bi[6][1] = -25.0/28.0 ,
453  bi[6][2] = 205.0/42.0 ,
454  bi[6][3] = -45.0/7.0 ,
455  bi[6][4] = 18.0/7.0 ,
456  // --------------------------------------------------------
457  //
458  // COEFFICIENTS OF bi[7]
459  bi[7][0] = 0.0 ,
460  bi[7][1] = -2.0/11.0 ,
461  bi[7][2] = 73.0/55.0 ,
462  bi[7][3] = -171.0/55.0 ,
463  bi[7][4] = 108.0/55.0 ,
464  // --------------------------------------------------------
465  //
466  // COEFFICIENTS OF bi[8]
467  bi[8][0] = 0.0 ,
468  bi[8][1] = 189.0/22.0 ,
469  bi[8][2] = -1593.0/55.0 ,
470  bi[8][3] = 3537.0/110.0 ,
471  bi[8][4] = -648.0/55.0 ,
472  // --------------------------------------------------------
473  //
474  // COEFFICIENTS OF bi[9]
475  bi[9][0] = 0.0 ,
476  bi[9][1] = 351.0/110.0 ,
477  bi[9][2] = -999.0/55.0 ,
478  bi[9][3] = 2943.0/110.0 ,
479  bi[9][4] = -648.0/55.0 ;
480  // --------------------------------------------------------
481 
482 
483 
484  for(G4int i = 0; i< numberOfVariables; i++)
485  yIn[i] = yInput[i];
486 
487  G4double tau0 = tau;
488  // Calculating the polynomials :
489 
490  for(int i=1; i<=9; i++){ //Here i is NOT the coordinate no. , it's stage no.
491  b[i] = 0;
492  tau = 1.0;
493  for(int j=0; j<=4; j++){
494  b[i] += bi[i][j]*tau;
495  tau*=tau0;
496  }
497  }
498 
499  for(int i=0; i<numberOfVariables; i++){ //Here i IS the cooridnate no.
500  yOut[i] = yIn[i] + Step*tau0*(b[1]*dydx[i] + b[2]*ak2[i] + b[3]*ak3[i] +
501  b[4]*ak4[i] + b[5]*ak5[i] + b[6]*ak6[i] +
502  b[7]*ak7[i] + b[8]*ak8[i] + b[9]*ak9[i] );
503  }
504 
505 }
506 
507 
508 
509 
510 
511 
512 
CLHEP::Hep3Vector G4ThreeVector
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
void RightHandSide(const double y[], double dydx[])
G4int GetNumberOfVariables() const
void SetupInterpolate(const G4double yInput[], const G4double dydx[], const G4double Step)
void Interpolate(const G4double yInput[], const G4double dydx[], const G4double Step, G4double yOut[], G4double tau)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4FSALDormandPrince745 * fAuxStepper
void interpolate(const G4double yInput[], const G4double dydx[], G4double yOut[], G4double Step, G4double tau)
G4FSALDormandPrince745(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
int G4int
Definition: G4Types.hh:78
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[], G4double nextDydx[])
Definition: Step.hh:41