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 | // |
---|
27 | // $Id: G4Integrator.icc,v 1.13 2006/06/29 18:59:47 gunter Exp $ |
---|
28 | // GEANT4 tag $Name: HEAD $ |
---|
29 | // |
---|
30 | // Implementation of G4Integrator methods. |
---|
31 | // |
---|
32 | // |
---|
33 | |
---|
34 | ///////////////////////////////////////////////////////////////////// |
---|
35 | // |
---|
36 | // Sympson integration method |
---|
37 | // |
---|
38 | ///////////////////////////////////////////////////////////////////// |
---|
39 | // |
---|
40 | // Integration of class member functions T::f by Simpson method. |
---|
41 | |
---|
42 | template <class T, class F> |
---|
43 | G4double G4Integrator<T,F>::Simpson( T& typeT, |
---|
44 | F f, |
---|
45 | G4double xInitial, |
---|
46 | G4double xFinal, |
---|
47 | G4int iterationNumber ) |
---|
48 | { |
---|
49 | G4int i ; |
---|
50 | G4double step = (xFinal - xInitial)/iterationNumber ; |
---|
51 | G4double x = xInitial ; |
---|
52 | G4double xPlus = xInitial + 0.5*step ; |
---|
53 | G4double mean = ( (typeT.*f)(xInitial) + (typeT.*f)(xFinal) )*0.5 ; |
---|
54 | G4double sum = (typeT.*f)(xPlus) ; |
---|
55 | |
---|
56 | for(i=1;i<iterationNumber;i++) |
---|
57 | { |
---|
58 | x += step ; |
---|
59 | xPlus += step ; |
---|
60 | mean += (typeT.*f)(x) ; |
---|
61 | sum += (typeT.*f)(xPlus) ; |
---|
62 | } |
---|
63 | mean += 2.0*sum ; |
---|
64 | |
---|
65 | return mean*step/3.0 ; |
---|
66 | } |
---|
67 | |
---|
68 | ///////////////////////////////////////////////////////////////////// |
---|
69 | // |
---|
70 | // Integration of class member functions T::f by Simpson method. |
---|
71 | // Convenient to use with 'this' pointer |
---|
72 | |
---|
73 | template <class T, class F> |
---|
74 | G4double G4Integrator<T,F>::Simpson( T* ptrT, |
---|
75 | F f, |
---|
76 | G4double xInitial, |
---|
77 | G4double xFinal, |
---|
78 | G4int iterationNumber ) |
---|
79 | { |
---|
80 | G4int i ; |
---|
81 | G4double step = (xFinal - xInitial)/iterationNumber ; |
---|
82 | G4double x = xInitial ; |
---|
83 | G4double xPlus = xInitial + 0.5*step ; |
---|
84 | G4double mean = ( (ptrT->*f)(xInitial) + (ptrT->*f)(xFinal) )*0.5 ; |
---|
85 | G4double sum = (ptrT->*f)(xPlus) ; |
---|
86 | |
---|
87 | for(i=1;i<iterationNumber;i++) |
---|
88 | { |
---|
89 | x += step ; |
---|
90 | xPlus += step ; |
---|
91 | mean += (ptrT->*f)(x) ; |
---|
92 | sum += (ptrT->*f)(xPlus) ; |
---|
93 | } |
---|
94 | mean += 2.0*sum ; |
---|
95 | |
---|
96 | return mean*step/3.0 ; |
---|
97 | } |
---|
98 | |
---|
99 | ///////////////////////////////////////////////////////////////////// |
---|
100 | // |
---|
101 | // Integration of class member functions T::f by Simpson method. |
---|
102 | // Convenient to use, when function f is defined in global scope, i.e. in main() |
---|
103 | // program |
---|
104 | |
---|
105 | template <class T, class F> |
---|
106 | G4double G4Integrator<T,F>::Simpson( G4double (*f)(G4double), |
---|
107 | G4double xInitial, |
---|
108 | G4double xFinal, |
---|
109 | G4int iterationNumber ) |
---|
110 | { |
---|
111 | G4int i ; |
---|
112 | G4double step = (xFinal - xInitial)/iterationNumber ; |
---|
113 | G4double x = xInitial ; |
---|
114 | G4double xPlus = xInitial + 0.5*step ; |
---|
115 | G4double mean = ( (*f)(xInitial) + (*f)(xFinal) )*0.5 ; |
---|
116 | G4double sum = (*f)(xPlus) ; |
---|
117 | |
---|
118 | for(i=1;i<iterationNumber;i++) |
---|
119 | { |
---|
120 | x += step ; |
---|
121 | xPlus += step ; |
---|
122 | mean += (*f)(x) ; |
---|
123 | sum += (*f)(xPlus) ; |
---|
124 | } |
---|
125 | mean += 2.0*sum ; |
---|
126 | |
---|
127 | return mean*step/3.0 ; |
---|
128 | } |
---|
129 | |
---|
130 | ////////////////////////////////////////////////////////////////////////// |
---|
131 | // |
---|
132 | // Adaptive Gauss method |
---|
133 | // |
---|
134 | ////////////////////////////////////////////////////////////////////////// |
---|
135 | // |
---|
136 | // |
---|
137 | |
---|
138 | template <class T, class F> |
---|
139 | G4double G4Integrator<T,F>::Gauss( T& typeT, F f, |
---|
140 | G4double xInitial, G4double xFinal ) |
---|
141 | { |
---|
142 | static G4double root = 1.0/std::sqrt(3.0) ; |
---|
143 | |
---|
144 | G4double xMean = (xInitial + xFinal)/2.0 ; |
---|
145 | G4double Step = (xFinal - xInitial)/2.0 ; |
---|
146 | G4double delta = Step*root ; |
---|
147 | G4double sum = ((typeT.*f)(xMean + delta) + |
---|
148 | (typeT.*f)(xMean - delta)) ; |
---|
149 | |
---|
150 | return sum*Step ; |
---|
151 | } |
---|
152 | |
---|
153 | ////////////////////////////////////////////////////////////////////// |
---|
154 | // |
---|
155 | // |
---|
156 | |
---|
157 | template <class T, class F> G4double |
---|
158 | G4Integrator<T,F>::Gauss( T* ptrT, F f, G4double a, G4double b ) |
---|
159 | { |
---|
160 | return Gauss(*ptrT,f,a,b) ; |
---|
161 | } |
---|
162 | |
---|
163 | /////////////////////////////////////////////////////////////////////// |
---|
164 | // |
---|
165 | // |
---|
166 | |
---|
167 | template <class T, class F> |
---|
168 | G4double G4Integrator<T,F>::Gauss( G4double (*f)(G4double), |
---|
169 | G4double xInitial, G4double xFinal) |
---|
170 | { |
---|
171 | static G4double root = 1.0/std::sqrt(3.0) ; |
---|
172 | |
---|
173 | G4double xMean = (xInitial + xFinal)/2.0 ; |
---|
174 | G4double Step = (xFinal - xInitial)/2.0 ; |
---|
175 | G4double delta = Step*root ; |
---|
176 | G4double sum = ( (*f)(xMean + delta) + (*f)(xMean - delta) ) ; |
---|
177 | |
---|
178 | return sum*Step ; |
---|
179 | } |
---|
180 | |
---|
181 | /////////////////////////////////////////////////////////////////////////// |
---|
182 | // |
---|
183 | // |
---|
184 | |
---|
185 | template <class T, class F> |
---|
186 | void G4Integrator<T,F>::AdaptGauss( T& typeT, F f, G4double xInitial, |
---|
187 | G4double xFinal, G4double fTolerance, |
---|
188 | G4double& sum, |
---|
189 | G4int& depth ) |
---|
190 | { |
---|
191 | if(depth > 100) |
---|
192 | { |
---|
193 | G4cout<<"G4Integrator<T,F>::AdaptGauss: WARNING !!!"<<G4endl ; |
---|
194 | G4cout<<"Function varies too rapidly to get stated accuracy in 100 steps " |
---|
195 | <<G4endl ; |
---|
196 | |
---|
197 | return ; |
---|
198 | } |
---|
199 | G4double xMean = (xInitial + xFinal)/2.0 ; |
---|
200 | G4double leftHalf = Gauss(typeT,f,xInitial,xMean) ; |
---|
201 | G4double rightHalf = Gauss(typeT,f,xMean,xFinal) ; |
---|
202 | G4double full = Gauss(typeT,f,xInitial,xFinal) ; |
---|
203 | if(std::fabs(leftHalf+rightHalf-full) < fTolerance) |
---|
204 | { |
---|
205 | sum += full ; |
---|
206 | } |
---|
207 | else |
---|
208 | { |
---|
209 | depth++ ; |
---|
210 | AdaptGauss(typeT,f,xInitial,xMean,fTolerance,sum,depth) ; |
---|
211 | AdaptGauss(typeT,f,xMean,xFinal,fTolerance,sum,depth) ; |
---|
212 | } |
---|
213 | } |
---|
214 | |
---|
215 | template <class T, class F> |
---|
216 | void G4Integrator<T,F>::AdaptGauss( T* ptrT, F f, G4double xInitial, |
---|
217 | G4double xFinal, G4double fTolerance, |
---|
218 | G4double& sum, |
---|
219 | G4int& depth ) |
---|
220 | { |
---|
221 | AdaptGauss(*ptrT,f,xInitial,xFinal,fTolerance,sum,depth) ; |
---|
222 | } |
---|
223 | |
---|
224 | ///////////////////////////////////////////////////////////////////////// |
---|
225 | // |
---|
226 | // |
---|
227 | template <class T, class F> |
---|
228 | void G4Integrator<T,F>::AdaptGauss( G4double (*f)(G4double), |
---|
229 | G4double xInitial, G4double xFinal, |
---|
230 | G4double fTolerance, G4double& sum, |
---|
231 | G4int& depth ) |
---|
232 | { |
---|
233 | if(depth > 100) |
---|
234 | { |
---|
235 | G4cout<<"G4SimpleIntegration::AdaptGauss: WARNING !!!"<<G4endl ; |
---|
236 | G4cout<<"Function varies too rapidly to get stated accuracy in 100 steps " |
---|
237 | <<G4endl ; |
---|
238 | |
---|
239 | return ; |
---|
240 | } |
---|
241 | G4double xMean = (xInitial + xFinal)/2.0 ; |
---|
242 | G4double leftHalf = Gauss(f,xInitial,xMean) ; |
---|
243 | G4double rightHalf = Gauss(f,xMean,xFinal) ; |
---|
244 | G4double full = Gauss(f,xInitial,xFinal) ; |
---|
245 | if(std::fabs(leftHalf+rightHalf-full) < fTolerance) |
---|
246 | { |
---|
247 | sum += full ; |
---|
248 | } |
---|
249 | else |
---|
250 | { |
---|
251 | depth++ ; |
---|
252 | AdaptGauss(f,xInitial,xMean,fTolerance,sum,depth) ; |
---|
253 | AdaptGauss(f,xMean,xFinal,fTolerance,sum,depth) ; |
---|
254 | } |
---|
255 | } |
---|
256 | |
---|
257 | //////////////////////////////////////////////////////////////////////// |
---|
258 | // |
---|
259 | // Adaptive Gauss integration with accuracy 'e' |
---|
260 | // Convenient for using with class object typeT |
---|
261 | |
---|
262 | template<class T, class F> |
---|
263 | G4double G4Integrator<T,F>::AdaptiveGauss( T& typeT, F f, G4double xInitial, |
---|
264 | G4double xFinal, G4double e ) |
---|
265 | { |
---|
266 | G4int depth = 0 ; |
---|
267 | G4double sum = 0.0 ; |
---|
268 | AdaptGauss(typeT,f,xInitial,xFinal,e,sum,depth) ; |
---|
269 | return sum ; |
---|
270 | } |
---|
271 | |
---|
272 | //////////////////////////////////////////////////////////////////////// |
---|
273 | // |
---|
274 | // Adaptive Gauss integration with accuracy 'e' |
---|
275 | // Convenient for using with 'this' pointer |
---|
276 | |
---|
277 | template<class T, class F> |
---|
278 | G4double G4Integrator<T,F>::AdaptiveGauss( T* ptrT, F f, G4double xInitial, |
---|
279 | G4double xFinal, G4double e ) |
---|
280 | { |
---|
281 | return AdaptiveGauss(*ptrT,f,xInitial,xFinal,e) ; |
---|
282 | } |
---|
283 | |
---|
284 | //////////////////////////////////////////////////////////////////////// |
---|
285 | // |
---|
286 | // Adaptive Gauss integration with accuracy 'e' |
---|
287 | // Convenient for using with global scope function f |
---|
288 | |
---|
289 | template <class T, class F> |
---|
290 | G4double G4Integrator<T,F>::AdaptiveGauss( G4double (*f)(G4double), |
---|
291 | G4double xInitial, G4double xFinal, G4double e ) |
---|
292 | { |
---|
293 | G4int depth = 0 ; |
---|
294 | G4double sum = 0.0 ; |
---|
295 | AdaptGauss(f,xInitial,xFinal,e,sum,depth) ; |
---|
296 | return sum ; |
---|
297 | } |
---|
298 | |
---|
299 | //////////////////////////////////////////////////////////////////////////// |
---|
300 | // Gauss integration methods involving ortogonal polynomials |
---|
301 | //////////////////////////////////////////////////////////////////////////// |
---|
302 | // |
---|
303 | // Methods involving Legendre polynomials |
---|
304 | // |
---|
305 | ///////////////////////////////////////////////////////////////////////// |
---|
306 | // |
---|
307 | // The value nLegendre set the accuracy required, i.e the number of points |
---|
308 | // where the function pFunction will be evaluated during integration. |
---|
309 | // The function creates the arrays for abscissas and weights that used |
---|
310 | // in Gauss-Legendre quadrature method. |
---|
311 | // The values a and b are the limits of integration of the function f . |
---|
312 | // nLegendre MUST BE EVEN !!! |
---|
313 | // Returns the integral of the function f between a and b, by 2*fNumber point |
---|
314 | // Gauss-Legendre integration: the function is evaluated exactly |
---|
315 | // 2*fNumber times at interior points in the range of integration. |
---|
316 | // Since the weights and abscissas are, in this case, symmetric around |
---|
317 | // the midpoint of the range of integration, there are actually only |
---|
318 | // fNumber distinct values of each. |
---|
319 | // Convenient for using with some class object dataT |
---|
320 | |
---|
321 | template <class T, class F> |
---|
322 | G4double G4Integrator<T,F>::Legendre( T& typeT, F f, G4double a, G4double b, |
---|
323 | G4int nLegendre ) |
---|
324 | { |
---|
325 | G4double newton, newton1, temp1, temp2, temp3, temp ; |
---|
326 | G4double xDiff, xMean, dx, integral ; |
---|
327 | |
---|
328 | const G4double tolerance = 1.6e-10 ; |
---|
329 | G4int i, j, k = nLegendre ; |
---|
330 | G4int fNumber = (nLegendre + 1)/2 ; |
---|
331 | |
---|
332 | if(2*fNumber != k) |
---|
333 | { |
---|
334 | G4Exception("G4Integrator<T,F>::Legendre(T&,F, ...)", "InvalidCall", |
---|
335 | FatalException, "Invalid (odd) nLegendre in constructor."); |
---|
336 | } |
---|
337 | |
---|
338 | G4double* fAbscissa = new G4double[fNumber] ; |
---|
339 | G4double* fWeight = new G4double[fNumber] ; |
---|
340 | |
---|
341 | for(i=1;i<=fNumber;i++) // Loop over the desired roots |
---|
342 | { |
---|
343 | newton = std::cos(pi*(i - 0.25)/(k + 0.5)) ; // Initial root approximation |
---|
344 | |
---|
345 | do // loop of Newton's method |
---|
346 | { |
---|
347 | temp1 = 1.0 ; |
---|
348 | temp2 = 0.0 ; |
---|
349 | for(j=1;j<=k;j++) |
---|
350 | { |
---|
351 | temp3 = temp2 ; |
---|
352 | temp2 = temp1 ; |
---|
353 | temp1 = ((2.0*j - 1.0)*newton*temp2 - (j - 1.0)*temp3)/j ; |
---|
354 | } |
---|
355 | temp = k*(newton*temp1 - temp2)/(newton*newton - 1.0) ; |
---|
356 | newton1 = newton ; |
---|
357 | newton = newton1 - temp1/temp ; // Newton's method |
---|
358 | } |
---|
359 | while(std::fabs(newton - newton1) > tolerance) ; |
---|
360 | |
---|
361 | fAbscissa[fNumber-i] = newton ; |
---|
362 | fWeight[fNumber-i] = 2.0/((1.0 - newton*newton)*temp*temp) ; |
---|
363 | } |
---|
364 | |
---|
365 | // |
---|
366 | // Now we ready to get integral |
---|
367 | // |
---|
368 | |
---|
369 | xMean = 0.5*(a + b) ; |
---|
370 | xDiff = 0.5*(b - a) ; |
---|
371 | integral = 0.0 ; |
---|
372 | for(i=0;i<fNumber;i++) |
---|
373 | { |
---|
374 | dx = xDiff*fAbscissa[i] ; |
---|
375 | integral += fWeight[i]*( (typeT.*f)(xMean + dx) + |
---|
376 | (typeT.*f)(xMean - dx) ) ; |
---|
377 | } |
---|
378 | delete[] fAbscissa; |
---|
379 | delete[] fWeight; |
---|
380 | return integral *= xDiff ; |
---|
381 | } |
---|
382 | |
---|
383 | /////////////////////////////////////////////////////////////////////// |
---|
384 | // |
---|
385 | // Convenient for using with the pointer 'this' |
---|
386 | |
---|
387 | template <class T, class F> |
---|
388 | G4double G4Integrator<T,F>::Legendre( T* ptrT, F f, G4double a, |
---|
389 | G4double b, G4int nLegendre ) |
---|
390 | { |
---|
391 | return Legendre(*ptrT,f,a,b,nLegendre) ; |
---|
392 | } |
---|
393 | |
---|
394 | /////////////////////////////////////////////////////////////////////// |
---|
395 | // |
---|
396 | // Convenient for using with global scope function f |
---|
397 | |
---|
398 | template <class T, class F> |
---|
399 | G4double G4Integrator<T,F>::Legendre( G4double (*f)(G4double), |
---|
400 | G4double a, G4double b, G4int nLegendre) |
---|
401 | { |
---|
402 | G4double newton, newton1, temp1, temp2, temp3, temp ; |
---|
403 | G4double xDiff, xMean, dx, integral ; |
---|
404 | |
---|
405 | const G4double tolerance = 1.6e-10 ; |
---|
406 | G4int i, j, k = nLegendre ; |
---|
407 | G4int fNumber = (nLegendre + 1)/2 ; |
---|
408 | |
---|
409 | if(2*fNumber != k) |
---|
410 | { |
---|
411 | G4Exception("G4Integrator<T,F>::Legendre(...)", "InvalidCall", |
---|
412 | FatalException, "Invalid (odd) nLegendre in constructor."); |
---|
413 | } |
---|
414 | |
---|
415 | G4double* fAbscissa = new G4double[fNumber] ; |
---|
416 | G4double* fWeight = new G4double[fNumber] ; |
---|
417 | |
---|
418 | for(i=1;i<=fNumber;i++) // Loop over the desired roots |
---|
419 | { |
---|
420 | newton = std::cos(pi*(i - 0.25)/(k + 0.5)) ; // Initial root approximation |
---|
421 | |
---|
422 | do // loop of Newton's method |
---|
423 | { |
---|
424 | temp1 = 1.0 ; |
---|
425 | temp2 = 0.0 ; |
---|
426 | for(j=1;j<=k;j++) |
---|
427 | { |
---|
428 | temp3 = temp2 ; |
---|
429 | temp2 = temp1 ; |
---|
430 | temp1 = ((2.0*j - 1.0)*newton*temp2 - (j - 1.0)*temp3)/j ; |
---|
431 | } |
---|
432 | temp = k*(newton*temp1 - temp2)/(newton*newton - 1.0) ; |
---|
433 | newton1 = newton ; |
---|
434 | newton = newton1 - temp1/temp ; // Newton's method |
---|
435 | } |
---|
436 | while(std::fabs(newton - newton1) > tolerance) ; |
---|
437 | |
---|
438 | fAbscissa[fNumber-i] = newton ; |
---|
439 | fWeight[fNumber-i] = 2.0/((1.0 - newton*newton)*temp*temp) ; |
---|
440 | } |
---|
441 | |
---|
442 | // |
---|
443 | // Now we ready to get integral |
---|
444 | // |
---|
445 | |
---|
446 | xMean = 0.5*(a + b) ; |
---|
447 | xDiff = 0.5*(b - a) ; |
---|
448 | integral = 0.0 ; |
---|
449 | for(i=0;i<fNumber;i++) |
---|
450 | { |
---|
451 | dx = xDiff*fAbscissa[i] ; |
---|
452 | integral += fWeight[i]*( (*f)(xMean + dx) + (*f)(xMean - dx) ) ; |
---|
453 | } |
---|
454 | delete[] fAbscissa; |
---|
455 | delete[] fWeight; |
---|
456 | |
---|
457 | return integral *= xDiff ; |
---|
458 | } |
---|
459 | |
---|
460 | //////////////////////////////////////////////////////////////////////////// |
---|
461 | // |
---|
462 | // Returns the integral of the function to be pointed by T::f between a and b, |
---|
463 | // by ten point Gauss-Legendre integration: the function is evaluated exactly |
---|
464 | // ten times at interior points in the range of integration. Since the weights |
---|
465 | // and abscissas are, in this case, symmetric around the midpoint of the |
---|
466 | // range of integration, there are actually only five distinct values of each |
---|
467 | // Convenient for using with class object typeT |
---|
468 | |
---|
469 | template <class T, class F> |
---|
470 | G4double G4Integrator<T,F>::Legendre10( T& typeT, F f,G4double a, G4double b) |
---|
471 | { |
---|
472 | G4int i ; |
---|
473 | G4double xDiff, xMean, dx, integral ; |
---|
474 | |
---|
475 | // From Abramowitz M., Stegan I.A. 1964 , Handbook of Math... , p. 916 |
---|
476 | |
---|
477 | static G4double abscissa[] = { 0.148874338981631, 0.433395394129247, |
---|
478 | 0.679409568299024, 0.865063366688985, |
---|
479 | 0.973906528517172 } ; |
---|
480 | |
---|
481 | static G4double weight[] = { 0.295524224714753, 0.269266719309996, |
---|
482 | 0.219086362515982, 0.149451349150581, |
---|
483 | 0.066671344308688 } ; |
---|
484 | xMean = 0.5*(a + b) ; |
---|
485 | xDiff = 0.5*(b - a) ; |
---|
486 | integral = 0.0 ; |
---|
487 | for(i=0;i<5;i++) |
---|
488 | { |
---|
489 | dx = xDiff*abscissa[i] ; |
---|
490 | integral += weight[i]*( (typeT.*f)(xMean + dx) + (typeT.*f)(xMean - dx)) ; |
---|
491 | } |
---|
492 | return integral *= xDiff ; |
---|
493 | } |
---|
494 | |
---|
495 | /////////////////////////////////////////////////////////////////////////// |
---|
496 | // |
---|
497 | // Convenient for using with the pointer 'this' |
---|
498 | |
---|
499 | template <class T, class F> |
---|
500 | G4double G4Integrator<T,F>::Legendre10( T* ptrT, F f,G4double a, G4double b) |
---|
501 | { |
---|
502 | return Legendre10(*ptrT,f,a,b) ; |
---|
503 | } |
---|
504 | |
---|
505 | ////////////////////////////////////////////////////////////////////////// |
---|
506 | // |
---|
507 | // Convenient for using with global scope functions |
---|
508 | |
---|
509 | template <class T, class F> |
---|
510 | G4double G4Integrator<T,F>::Legendre10( G4double (*f)(G4double), |
---|
511 | G4double a, G4double b ) |
---|
512 | { |
---|
513 | G4int i ; |
---|
514 | G4double xDiff, xMean, dx, integral ; |
---|
515 | |
---|
516 | // From Abramowitz M., Stegan I.A. 1964 , Handbook of Math... , p. 916 |
---|
517 | |
---|
518 | static G4double abscissa[] = { 0.148874338981631, 0.433395394129247, |
---|
519 | 0.679409568299024, 0.865063366688985, |
---|
520 | 0.973906528517172 } ; |
---|
521 | |
---|
522 | static G4double weight[] = { 0.295524224714753, 0.269266719309996, |
---|
523 | 0.219086362515982, 0.149451349150581, |
---|
524 | 0.066671344308688 } ; |
---|
525 | xMean = 0.5*(a + b) ; |
---|
526 | xDiff = 0.5*(b - a) ; |
---|
527 | integral = 0.0 ; |
---|
528 | for(i=0;i<5;i++) |
---|
529 | { |
---|
530 | dx = xDiff*abscissa[i] ; |
---|
531 | integral += weight[i]*( (*f)(xMean + dx) + (*f)(xMean - dx)) ; |
---|
532 | } |
---|
533 | return integral *= xDiff ; |
---|
534 | } |
---|
535 | |
---|
536 | /////////////////////////////////////////////////////////////////////// |
---|
537 | // |
---|
538 | // Returns the integral of the function to be pointed by T::f between a and b, |
---|
539 | // by 96 point Gauss-Legendre integration: the function is evaluated exactly |
---|
540 | // ten Times at interior points in the range of integration. Since the weights |
---|
541 | // and abscissas are, in this case, symmetric around the midpoint of the |
---|
542 | // range of integration, there are actually only five distinct values of each |
---|
543 | // Convenient for using with some class object typeT |
---|
544 | |
---|
545 | template <class T, class F> |
---|
546 | G4double G4Integrator<T,F>::Legendre96( T& typeT, F f,G4double a, G4double b) |
---|
547 | { |
---|
548 | G4int i ; |
---|
549 | G4double xDiff, xMean, dx, integral ; |
---|
550 | |
---|
551 | // From Abramowitz M., Stegan I.A. 1964 , Handbook of Math... , p. 919 |
---|
552 | |
---|
553 | static G4double |
---|
554 | abscissa[] = { |
---|
555 | 0.016276744849602969579, 0.048812985136049731112, |
---|
556 | 0.081297495464425558994, 0.113695850110665920911, |
---|
557 | 0.145973714654896941989, 0.178096882367618602759, // 6 |
---|
558 | |
---|
559 | 0.210031310460567203603, 0.241743156163840012328, |
---|
560 | 0.273198812591049141487, 0.304364944354496353024, |
---|
561 | 0.335208522892625422616, 0.365696861472313635031, // 12 |
---|
562 | |
---|
563 | 0.395797649828908603285, 0.425478988407300545365, |
---|
564 | 0.454709422167743008636, 0.483457973920596359768, |
---|
565 | 0.511694177154667673586, 0.539388108324357436227, // 18 |
---|
566 | |
---|
567 | 0.566510418561397168404, 0.593032364777572080684, |
---|
568 | 0.618925840125468570386, 0.644163403784967106798, |
---|
569 | 0.668718310043916153953, 0.692564536642171561344, // 24 |
---|
570 | |
---|
571 | 0.715676812348967626225, 0.738030643744400132851, |
---|
572 | 0.759602341176647498703, 0.780369043867433217604, |
---|
573 | 0.800308744139140817229, 0.819400310737931675539, // 30 |
---|
574 | |
---|
575 | 0.837623511228187121494, 0.854959033434601455463, |
---|
576 | 0.871388505909296502874, 0.886894517402420416057, |
---|
577 | 0.901460635315852341319, 0.915071423120898074206, // 36 |
---|
578 | |
---|
579 | 0.927712456722308690965, 0.939370339752755216932, |
---|
580 | 0.950032717784437635756, 0.959688291448742539300, |
---|
581 | 0.968326828463264212174, 0.975939174585136466453, // 42 |
---|
582 | |
---|
583 | 0.982517263563014677447, 0.988054126329623799481, |
---|
584 | 0.992543900323762624572, 0.995981842987209290650, |
---|
585 | 0.998364375863181677724, 0.999689503883230766828 // 48 |
---|
586 | } ; |
---|
587 | |
---|
588 | static G4double |
---|
589 | weight[] = { |
---|
590 | 0.032550614492363166242, 0.032516118713868835987, |
---|
591 | 0.032447163714064269364, 0.032343822568575928429, |
---|
592 | 0.032206204794030250669, 0.032034456231992663218, // 6 |
---|
593 | |
---|
594 | 0.031828758894411006535, 0.031589330770727168558, |
---|
595 | 0.031316425596862355813, 0.031010332586313837423, |
---|
596 | 0.030671376123669149014, 0.030299915420827593794, // 12 |
---|
597 | |
---|
598 | 0.029896344136328385984, 0.029461089958167905970, |
---|
599 | 0.028994614150555236543, 0.028497411065085385646, |
---|
600 | 0.027970007616848334440, 0.027412962726029242823, // 18 |
---|
601 | |
---|
602 | 0.026826866725591762198, 0.026212340735672413913, |
---|
603 | 0.025570036005349361499, 0.024900633222483610288, |
---|
604 | 0.024204841792364691282, 0.023483399085926219842, // 24 |
---|
605 | |
---|
606 | 0.022737069658329374001, 0.021966644438744349195, |
---|
607 | 0.021172939892191298988, 0.020356797154333324595, |
---|
608 | 0.019519081140145022410, 0.018660679627411467385, // 30 |
---|
609 | |
---|
610 | 0.017782502316045260838, 0.016885479864245172450, |
---|
611 | 0.015970562902562291381, 0.015038721026994938006, |
---|
612 | 0.014090941772314860916, 0.013128229566961572637, // 36 |
---|
613 | |
---|
614 | 0.012151604671088319635, 0.011162102099838498591, |
---|
615 | 0.010160770535008415758, 0.009148671230783386633, |
---|
616 | 0.008126876925698759217, 0.007096470791153865269, // 42 |
---|
617 | |
---|
618 | 0.006058545504235961683, 0.005014202742927517693, |
---|
619 | 0.003964554338444686674, 0.002910731817934946408, |
---|
620 | 0.001853960788946921732, 0.000796792065552012429 // 48 |
---|
621 | } ; |
---|
622 | xMean = 0.5*(a + b) ; |
---|
623 | xDiff = 0.5*(b - a) ; |
---|
624 | integral = 0.0 ; |
---|
625 | for(i=0;i<48;i++) |
---|
626 | { |
---|
627 | dx = xDiff*abscissa[i] ; |
---|
628 | integral += weight[i]*((typeT.*f)(xMean + dx) + (typeT.*f)(xMean - dx)) ; |
---|
629 | } |
---|
630 | return integral *= xDiff ; |
---|
631 | } |
---|
632 | |
---|
633 | /////////////////////////////////////////////////////////////////////// |
---|
634 | // |
---|
635 | // Convenient for using with the pointer 'this' |
---|
636 | |
---|
637 | template <class T, class F> |
---|
638 | G4double G4Integrator<T,F>::Legendre96( T* ptrT, F f,G4double a, G4double b) |
---|
639 | { |
---|
640 | return Legendre96(*ptrT,f,a,b) ; |
---|
641 | } |
---|
642 | |
---|
643 | /////////////////////////////////////////////////////////////////////// |
---|
644 | // |
---|
645 | // Convenient for using with global scope function f |
---|
646 | |
---|
647 | template <class T, class F> |
---|
648 | G4double G4Integrator<T,F>::Legendre96( G4double (*f)(G4double), |
---|
649 | G4double a, G4double b ) |
---|
650 | { |
---|
651 | G4int i ; |
---|
652 | G4double xDiff, xMean, dx, integral ; |
---|
653 | |
---|
654 | // From Abramowitz M., Stegan I.A. 1964 , Handbook of Math... , p. 919 |
---|
655 | |
---|
656 | static G4double |
---|
657 | abscissa[] = { |
---|
658 | 0.016276744849602969579, 0.048812985136049731112, |
---|
659 | 0.081297495464425558994, 0.113695850110665920911, |
---|
660 | 0.145973714654896941989, 0.178096882367618602759, // 6 |
---|
661 | |
---|
662 | 0.210031310460567203603, 0.241743156163840012328, |
---|
663 | 0.273198812591049141487, 0.304364944354496353024, |
---|
664 | 0.335208522892625422616, 0.365696861472313635031, // 12 |
---|
665 | |
---|
666 | 0.395797649828908603285, 0.425478988407300545365, |
---|
667 | 0.454709422167743008636, 0.483457973920596359768, |
---|
668 | 0.511694177154667673586, 0.539388108324357436227, // 18 |
---|
669 | |
---|
670 | 0.566510418561397168404, 0.593032364777572080684, |
---|
671 | 0.618925840125468570386, 0.644163403784967106798, |
---|
672 | 0.668718310043916153953, 0.692564536642171561344, // 24 |
---|
673 | |
---|
674 | 0.715676812348967626225, 0.738030643744400132851, |
---|
675 | 0.759602341176647498703, 0.780369043867433217604, |
---|
676 | 0.800308744139140817229, 0.819400310737931675539, // 30 |
---|
677 | |
---|
678 | 0.837623511228187121494, 0.854959033434601455463, |
---|
679 | 0.871388505909296502874, 0.886894517402420416057, |
---|
680 | 0.901460635315852341319, 0.915071423120898074206, // 36 |
---|
681 | |
---|
682 | 0.927712456722308690965, 0.939370339752755216932, |
---|
683 | 0.950032717784437635756, 0.959688291448742539300, |
---|
684 | 0.968326828463264212174, 0.975939174585136466453, // 42 |
---|
685 | |
---|
686 | 0.982517263563014677447, 0.988054126329623799481, |
---|
687 | 0.992543900323762624572, 0.995981842987209290650, |
---|
688 | 0.998364375863181677724, 0.999689503883230766828 // 48 |
---|
689 | } ; |
---|
690 | |
---|
691 | static G4double |
---|
692 | weight[] = { |
---|
693 | 0.032550614492363166242, 0.032516118713868835987, |
---|
694 | 0.032447163714064269364, 0.032343822568575928429, |
---|
695 | 0.032206204794030250669, 0.032034456231992663218, // 6 |
---|
696 | |
---|
697 | 0.031828758894411006535, 0.031589330770727168558, |
---|
698 | 0.031316425596862355813, 0.031010332586313837423, |
---|
699 | 0.030671376123669149014, 0.030299915420827593794, // 12 |
---|
700 | |
---|
701 | 0.029896344136328385984, 0.029461089958167905970, |
---|
702 | 0.028994614150555236543, 0.028497411065085385646, |
---|
703 | 0.027970007616848334440, 0.027412962726029242823, // 18 |
---|
704 | |
---|
705 | 0.026826866725591762198, 0.026212340735672413913, |
---|
706 | 0.025570036005349361499, 0.024900633222483610288, |
---|
707 | 0.024204841792364691282, 0.023483399085926219842, // 24 |
---|
708 | |
---|
709 | 0.022737069658329374001, 0.021966644438744349195, |
---|
710 | 0.021172939892191298988, 0.020356797154333324595, |
---|
711 | 0.019519081140145022410, 0.018660679627411467385, // 30 |
---|
712 | |
---|
713 | 0.017782502316045260838, 0.016885479864245172450, |
---|
714 | 0.015970562902562291381, 0.015038721026994938006, |
---|
715 | 0.014090941772314860916, 0.013128229566961572637, // 36 |
---|
716 | |
---|
717 | 0.012151604671088319635, 0.011162102099838498591, |
---|
718 | 0.010160770535008415758, 0.009148671230783386633, |
---|
719 | 0.008126876925698759217, 0.007096470791153865269, // 42 |
---|
720 | |
---|
721 | 0.006058545504235961683, 0.005014202742927517693, |
---|
722 | 0.003964554338444686674, 0.002910731817934946408, |
---|
723 | 0.001853960788946921732, 0.000796792065552012429 // 48 |
---|
724 | } ; |
---|
725 | xMean = 0.5*(a + b) ; |
---|
726 | xDiff = 0.5*(b - a) ; |
---|
727 | integral = 0.0 ; |
---|
728 | for(i=0;i<48;i++) |
---|
729 | { |
---|
730 | dx = xDiff*abscissa[i] ; |
---|
731 | integral += weight[i]*((*f)(xMean + dx) + (*f)(xMean - dx)) ; |
---|
732 | } |
---|
733 | return integral *= xDiff ; |
---|
734 | } |
---|
735 | |
---|
736 | ////////////////////////////////////////////////////////////////////////////// |
---|
737 | // |
---|
738 | // Methods involving Chebyshev polynomials |
---|
739 | // |
---|
740 | /////////////////////////////////////////////////////////////////////////// |
---|
741 | // |
---|
742 | // Integrates function pointed by T::f from a to b by Gauss-Chebyshev |
---|
743 | // quadrature method. |
---|
744 | // Convenient for using with class object typeT |
---|
745 | |
---|
746 | template <class T, class F> |
---|
747 | G4double G4Integrator<T,F>::Chebyshev( T& typeT, F f, G4double a, |
---|
748 | G4double b, G4int nChebyshev ) |
---|
749 | { |
---|
750 | G4int i ; |
---|
751 | G4double xDiff, xMean, dx, integral = 0.0 ; |
---|
752 | |
---|
753 | G4int fNumber = nChebyshev ; // Try to reduce fNumber twice ?? |
---|
754 | G4double cof = pi/fNumber ; |
---|
755 | G4double* fAbscissa = new G4double[fNumber] ; |
---|
756 | G4double* fWeight = new G4double[fNumber] ; |
---|
757 | for(i=0;i<fNumber;i++) |
---|
758 | { |
---|
759 | fAbscissa[i] = std::cos(cof*(i + 0.5)) ; |
---|
760 | fWeight[i] = cof*std::sqrt(1 - fAbscissa[i]*fAbscissa[i]) ; |
---|
761 | } |
---|
762 | |
---|
763 | // |
---|
764 | // Now we ready to estimate the integral |
---|
765 | // |
---|
766 | |
---|
767 | xMean = 0.5*(a + b) ; |
---|
768 | xDiff = 0.5*(b - a) ; |
---|
769 | for(i=0;i<fNumber;i++) |
---|
770 | { |
---|
771 | dx = xDiff*fAbscissa[i] ; |
---|
772 | integral += fWeight[i]*(typeT.*f)(xMean + dx) ; |
---|
773 | } |
---|
774 | delete[] fAbscissa; |
---|
775 | delete[] fWeight; |
---|
776 | return integral *= xDiff ; |
---|
777 | } |
---|
778 | |
---|
779 | /////////////////////////////////////////////////////////////////////// |
---|
780 | // |
---|
781 | // Convenient for using with 'this' pointer |
---|
782 | |
---|
783 | template <class T, class F> |
---|
784 | G4double G4Integrator<T,F>::Chebyshev( T* ptrT, F f, G4double a, |
---|
785 | G4double b, G4int n ) |
---|
786 | { |
---|
787 | return Chebyshev(*ptrT,f,a,b,n) ; |
---|
788 | } |
---|
789 | |
---|
790 | //////////////////////////////////////////////////////////////////////// |
---|
791 | // |
---|
792 | // For use with global scope functions f |
---|
793 | |
---|
794 | template <class T, class F> |
---|
795 | G4double G4Integrator<T,F>::Chebyshev( G4double (*f)(G4double), |
---|
796 | G4double a, G4double b, G4int nChebyshev ) |
---|
797 | { |
---|
798 | G4int i ; |
---|
799 | G4double xDiff, xMean, dx, integral = 0.0 ; |
---|
800 | |
---|
801 | G4int fNumber = nChebyshev ; // Try to reduce fNumber twice ?? |
---|
802 | G4double cof = pi/fNumber ; |
---|
803 | G4double* fAbscissa = new G4double[fNumber] ; |
---|
804 | G4double* fWeight = new G4double[fNumber] ; |
---|
805 | for(i=0;i<fNumber;i++) |
---|
806 | { |
---|
807 | fAbscissa[i] = std::cos(cof*(i + 0.5)) ; |
---|
808 | fWeight[i] = cof*std::sqrt(1 - fAbscissa[i]*fAbscissa[i]) ; |
---|
809 | } |
---|
810 | |
---|
811 | // |
---|
812 | // Now we ready to estimate the integral |
---|
813 | // |
---|
814 | |
---|
815 | xMean = 0.5*(a + b) ; |
---|
816 | xDiff = 0.5*(b - a) ; |
---|
817 | for(i=0;i<fNumber;i++) |
---|
818 | { |
---|
819 | dx = xDiff*fAbscissa[i] ; |
---|
820 | integral += fWeight[i]*(*f)(xMean + dx) ; |
---|
821 | } |
---|
822 | delete[] fAbscissa; |
---|
823 | delete[] fWeight; |
---|
824 | return integral *= xDiff ; |
---|
825 | } |
---|
826 | |
---|
827 | ////////////////////////////////////////////////////////////////////// |
---|
828 | // |
---|
829 | // Method involving Laguerre polynomials |
---|
830 | // |
---|
831 | ////////////////////////////////////////////////////////////////////// |
---|
832 | // |
---|
833 | // Integral from zero to infinity of std::pow(x,alpha)*std::exp(-x)*f(x). |
---|
834 | // The value of nLaguerre sets the accuracy. |
---|
835 | // The function creates arrays fAbscissa[0,..,nLaguerre-1] and |
---|
836 | // fWeight[0,..,nLaguerre-1] . |
---|
837 | // Convenient for using with class object 'typeT' and (typeT.*f) function |
---|
838 | // (T::f) |
---|
839 | |
---|
840 | template <class T, class F> |
---|
841 | G4double G4Integrator<T,F>::Laguerre( T& typeT, F f, G4double alpha, |
---|
842 | G4int nLaguerre ) |
---|
843 | { |
---|
844 | const G4double tolerance = 1.0e-10 ; |
---|
845 | const G4int maxNumber = 12 ; |
---|
846 | G4int i, j, k ; |
---|
847 | G4double newton=0., newton1, temp1, temp2, temp3, temp, cofi ; |
---|
848 | G4double integral = 0.0 ; |
---|
849 | |
---|
850 | G4int fNumber = nLaguerre ; |
---|
851 | G4double* fAbscissa = new G4double[fNumber] ; |
---|
852 | G4double* fWeight = new G4double[fNumber] ; |
---|
853 | |
---|
854 | for(i=1;i<=fNumber;i++) // Loop over the desired roots |
---|
855 | { |
---|
856 | if(i == 1) |
---|
857 | { |
---|
858 | newton = (1.0 + alpha)*(3.0 + 0.92*alpha) |
---|
859 | / (1.0 + 2.4*fNumber + 1.8*alpha) ; |
---|
860 | } |
---|
861 | else if(i == 2) |
---|
862 | { |
---|
863 | newton += (15.0 + 6.25*alpha)/(1.0 + 0.9*alpha + 2.5*fNumber) ; |
---|
864 | } |
---|
865 | else |
---|
866 | { |
---|
867 | cofi = i - 2 ; |
---|
868 | newton += ((1.0+2.55*cofi)/(1.9*cofi) |
---|
869 | + 1.26*cofi*alpha/(1.0+3.5*cofi)) |
---|
870 | * (newton - fAbscissa[i-3])/(1.0 + 0.3*alpha) ; |
---|
871 | } |
---|
872 | for(k=1;k<=maxNumber;k++) |
---|
873 | { |
---|
874 | temp1 = 1.0 ; |
---|
875 | temp2 = 0.0 ; |
---|
876 | |
---|
877 | for(j=1;j<=fNumber;j++) |
---|
878 | { |
---|
879 | temp3 = temp2 ; |
---|
880 | temp2 = temp1 ; |
---|
881 | temp1 = ((2*j - 1 + alpha - newton)*temp2 - (j - 1 + alpha)*temp3)/j ; |
---|
882 | } |
---|
883 | temp = (fNumber*temp1 - (fNumber +alpha)*temp2)/newton ; |
---|
884 | newton1 = newton ; |
---|
885 | newton = newton1 - temp1/temp ; |
---|
886 | |
---|
887 | if(std::fabs(newton - newton1) <= tolerance) |
---|
888 | { |
---|
889 | break ; |
---|
890 | } |
---|
891 | } |
---|
892 | if(k > maxNumber) |
---|
893 | { |
---|
894 | G4Exception("G4Integrator<T,F>::Laguerre(T,F, ...)", "Error", |
---|
895 | FatalException, "Too many (>12) iterations."); |
---|
896 | } |
---|
897 | |
---|
898 | fAbscissa[i-1] = newton ; |
---|
899 | fWeight[i-1] = -std::exp(GammaLogarithm(alpha + fNumber) - |
---|
900 | GammaLogarithm((G4double)fNumber))/(temp*fNumber*temp2) ; |
---|
901 | } |
---|
902 | |
---|
903 | // |
---|
904 | // Integral evaluation |
---|
905 | // |
---|
906 | |
---|
907 | for(i=0;i<fNumber;i++) |
---|
908 | { |
---|
909 | integral += fWeight[i]*(typeT.*f)(fAbscissa[i]) ; |
---|
910 | } |
---|
911 | delete[] fAbscissa; |
---|
912 | delete[] fWeight; |
---|
913 | return integral ; |
---|
914 | } |
---|
915 | |
---|
916 | |
---|
917 | |
---|
918 | ////////////////////////////////////////////////////////////////////// |
---|
919 | // |
---|
920 | // |
---|
921 | |
---|
922 | template <class T, class F> G4double |
---|
923 | G4Integrator<T,F>::Laguerre( T* ptrT, F f, G4double alpha, G4int nLaguerre ) |
---|
924 | { |
---|
925 | return Laguerre(*ptrT,f,alpha,nLaguerre) ; |
---|
926 | } |
---|
927 | |
---|
928 | //////////////////////////////////////////////////////////////////////// |
---|
929 | // |
---|
930 | // For use with global scope functions f |
---|
931 | |
---|
932 | template <class T, class F> G4double |
---|
933 | G4Integrator<T,F>::Laguerre( G4double (*f)(G4double), |
---|
934 | G4double alpha, G4int nLaguerre ) |
---|
935 | { |
---|
936 | const G4double tolerance = 1.0e-10 ; |
---|
937 | const G4int maxNumber = 12 ; |
---|
938 | G4int i, j, k ; |
---|
939 | G4double newton=0., newton1, temp1, temp2, temp3, temp, cofi ; |
---|
940 | G4double integral = 0.0 ; |
---|
941 | |
---|
942 | G4int fNumber = nLaguerre ; |
---|
943 | G4double* fAbscissa = new G4double[fNumber] ; |
---|
944 | G4double* fWeight = new G4double[fNumber] ; |
---|
945 | |
---|
946 | for(i=1;i<=fNumber;i++) // Loop over the desired roots |
---|
947 | { |
---|
948 | if(i == 1) |
---|
949 | { |
---|
950 | newton = (1.0 + alpha)*(3.0 + 0.92*alpha) |
---|
951 | / (1.0 + 2.4*fNumber + 1.8*alpha) ; |
---|
952 | } |
---|
953 | else if(i == 2) |
---|
954 | { |
---|
955 | newton += (15.0 + 6.25*alpha)/(1.0 + 0.9*alpha + 2.5*fNumber) ; |
---|
956 | } |
---|
957 | else |
---|
958 | { |
---|
959 | cofi = i - 2 ; |
---|
960 | newton += ((1.0+2.55*cofi)/(1.9*cofi) |
---|
961 | + 1.26*cofi*alpha/(1.0+3.5*cofi)) |
---|
962 | * (newton - fAbscissa[i-3])/(1.0 + 0.3*alpha) ; |
---|
963 | } |
---|
964 | for(k=1;k<=maxNumber;k++) |
---|
965 | { |
---|
966 | temp1 = 1.0 ; |
---|
967 | temp2 = 0.0 ; |
---|
968 | |
---|
969 | for(j=1;j<=fNumber;j++) |
---|
970 | { |
---|
971 | temp3 = temp2 ; |
---|
972 | temp2 = temp1 ; |
---|
973 | temp1 = ((2*j - 1 + alpha - newton)*temp2 - (j - 1 + alpha)*temp3)/j ; |
---|
974 | } |
---|
975 | temp = (fNumber*temp1 - (fNumber +alpha)*temp2)/newton ; |
---|
976 | newton1 = newton ; |
---|
977 | newton = newton1 - temp1/temp ; |
---|
978 | |
---|
979 | if(std::fabs(newton - newton1) <= tolerance) |
---|
980 | { |
---|
981 | break ; |
---|
982 | } |
---|
983 | } |
---|
984 | if(k > maxNumber) |
---|
985 | { |
---|
986 | G4Exception("G4Integrator<T,F>::Laguerre( ...)", "Error", |
---|
987 | FatalException, "Too many (>12) iterations."); |
---|
988 | } |
---|
989 | |
---|
990 | fAbscissa[i-1] = newton ; |
---|
991 | fWeight[i-1] = -std::exp(GammaLogarithm(alpha + fNumber) - |
---|
992 | GammaLogarithm((G4double)fNumber))/(temp*fNumber*temp2) ; |
---|
993 | } |
---|
994 | |
---|
995 | // |
---|
996 | // Integral evaluation |
---|
997 | // |
---|
998 | |
---|
999 | for(i=0;i<fNumber;i++) |
---|
1000 | { |
---|
1001 | integral += fWeight[i]*(*f)(fAbscissa[i]) ; |
---|
1002 | } |
---|
1003 | delete[] fAbscissa; |
---|
1004 | delete[] fWeight; |
---|
1005 | return integral ; |
---|
1006 | } |
---|
1007 | |
---|
1008 | /////////////////////////////////////////////////////////////////////// |
---|
1009 | // |
---|
1010 | // Auxiliary function which returns the value of std::log(gamma-function(x)) |
---|
1011 | // Returns the value ln(Gamma(xx) for xx > 0. Full accuracy is obtained for |
---|
1012 | // xx > 1. For 0 < xx < 1. the reflection formula (6.1.4) can be used first. |
---|
1013 | // (Adapted from Numerical Recipes in C) |
---|
1014 | // |
---|
1015 | |
---|
1016 | template <class T, class F> |
---|
1017 | G4double G4Integrator<T,F>::GammaLogarithm(G4double xx) |
---|
1018 | { |
---|
1019 | static G4double cof[6] = { 76.18009172947146, -86.50532032941677, |
---|
1020 | 24.01409824083091, -1.231739572450155, |
---|
1021 | 0.1208650973866179e-2, -0.5395239384953e-5 } ; |
---|
1022 | register G4int j; |
---|
1023 | G4double x = xx - 1.0 ; |
---|
1024 | G4double tmp = x + 5.5 ; |
---|
1025 | tmp -= (x + 0.5) * std::log(tmp) ; |
---|
1026 | G4double ser = 1.000000000190015 ; |
---|
1027 | |
---|
1028 | for ( j = 0; j <= 5; j++ ) |
---|
1029 | { |
---|
1030 | x += 1.0 ; |
---|
1031 | ser += cof[j]/x ; |
---|
1032 | } |
---|
1033 | return -tmp + std::log(2.5066282746310005*ser) ; |
---|
1034 | } |
---|
1035 | |
---|
1036 | /////////////////////////////////////////////////////////////////////// |
---|
1037 | // |
---|
1038 | // Method involving Hermite polynomials |
---|
1039 | // |
---|
1040 | /////////////////////////////////////////////////////////////////////// |
---|
1041 | // |
---|
1042 | // |
---|
1043 | // Gauss-Hermite method for integration of std::exp(-x*x)*f(x) |
---|
1044 | // from minus infinity to plus infinity . |
---|
1045 | // |
---|
1046 | |
---|
1047 | template <class T, class F> |
---|
1048 | G4double G4Integrator<T,F>::Hermite( T& typeT, F f, G4int nHermite ) |
---|
1049 | { |
---|
1050 | const G4double tolerance = 1.0e-12 ; |
---|
1051 | const G4int maxNumber = 12 ; |
---|
1052 | |
---|
1053 | G4int i, j, k ; |
---|
1054 | G4double integral = 0.0 ; |
---|
1055 | G4double newton=0., newton1, temp1, temp2, temp3, temp ; |
---|
1056 | |
---|
1057 | G4double piInMinusQ = std::pow(pi,-0.25) ; // 1.0/std::sqrt(std::sqrt(pi)) ?? |
---|
1058 | |
---|
1059 | G4int fNumber = (nHermite +1)/2 ; |
---|
1060 | G4double* fAbscissa = new G4double[fNumber] ; |
---|
1061 | G4double* fWeight = new G4double[fNumber] ; |
---|
1062 | |
---|
1063 | for(i=1;i<=fNumber;i++) |
---|
1064 | { |
---|
1065 | if(i == 1) |
---|
1066 | { |
---|
1067 | newton = std::sqrt((G4double)(2*nHermite + 1)) - |
---|
1068 | 1.85575001*std::pow((G4double)(2*nHermite + 1),-0.16666999) ; |
---|
1069 | } |
---|
1070 | else if(i == 2) |
---|
1071 | { |
---|
1072 | newton -= 1.14001*std::pow((G4double)nHermite,0.425999)/newton ; |
---|
1073 | } |
---|
1074 | else if(i == 3) |
---|
1075 | { |
---|
1076 | newton = 1.86002*newton - 0.86002*fAbscissa[0] ; |
---|
1077 | } |
---|
1078 | else if(i == 4) |
---|
1079 | { |
---|
1080 | newton = 1.91001*newton - 0.91001*fAbscissa[1] ; |
---|
1081 | } |
---|
1082 | else |
---|
1083 | { |
---|
1084 | newton = 2.0*newton - fAbscissa[i - 3] ; |
---|
1085 | } |
---|
1086 | for(k=1;k<=maxNumber;k++) |
---|
1087 | { |
---|
1088 | temp1 = piInMinusQ ; |
---|
1089 | temp2 = 0.0 ; |
---|
1090 | |
---|
1091 | for(j=1;j<=nHermite;j++) |
---|
1092 | { |
---|
1093 | temp3 = temp2 ; |
---|
1094 | temp2 = temp1 ; |
---|
1095 | temp1 = newton*std::sqrt(2.0/j)*temp2 - |
---|
1096 | std::sqrt(((G4double)(j - 1))/j)*temp3 ; |
---|
1097 | } |
---|
1098 | temp = std::sqrt((G4double)2*nHermite)*temp2 ; |
---|
1099 | newton1 = newton ; |
---|
1100 | newton = newton1 - temp1/temp ; |
---|
1101 | |
---|
1102 | if(std::fabs(newton - newton1) <= tolerance) |
---|
1103 | { |
---|
1104 | break ; |
---|
1105 | } |
---|
1106 | } |
---|
1107 | if(k > maxNumber) |
---|
1108 | { |
---|
1109 | G4Exception("G4Integrator<T,F>::Hermite(T,F, ...)", "Error", |
---|
1110 | FatalException, "Too many (>12) iterations."); |
---|
1111 | } |
---|
1112 | fAbscissa[i-1] = newton ; |
---|
1113 | fWeight[i-1] = 2.0/(temp*temp) ; |
---|
1114 | } |
---|
1115 | |
---|
1116 | // |
---|
1117 | // Integral calculation |
---|
1118 | // |
---|
1119 | |
---|
1120 | for(i=0;i<fNumber;i++) |
---|
1121 | { |
---|
1122 | integral += fWeight[i]*( (typeT.*f)(fAbscissa[i]) + |
---|
1123 | (typeT.*f)(-fAbscissa[i]) ) ; |
---|
1124 | } |
---|
1125 | delete[] fAbscissa; |
---|
1126 | delete[] fWeight; |
---|
1127 | return integral ; |
---|
1128 | } |
---|
1129 | |
---|
1130 | |
---|
1131 | //////////////////////////////////////////////////////////////////////// |
---|
1132 | // |
---|
1133 | // For use with 'this' pointer |
---|
1134 | |
---|
1135 | template <class T, class F> |
---|
1136 | G4double G4Integrator<T,F>::Hermite( T* ptrT, F f, G4int n ) |
---|
1137 | { |
---|
1138 | return Hermite(*ptrT,f,n) ; |
---|
1139 | } |
---|
1140 | |
---|
1141 | //////////////////////////////////////////////////////////////////////// |
---|
1142 | // |
---|
1143 | // For use with global scope f |
---|
1144 | |
---|
1145 | template <class T, class F> |
---|
1146 | G4double G4Integrator<T,F>::Hermite( G4double (*f)(G4double), G4int nHermite) |
---|
1147 | { |
---|
1148 | const G4double tolerance = 1.0e-12 ; |
---|
1149 | const G4int maxNumber = 12 ; |
---|
1150 | |
---|
1151 | G4int i, j, k ; |
---|
1152 | G4double integral = 0.0 ; |
---|
1153 | G4double newton=0., newton1, temp1, temp2, temp3, temp ; |
---|
1154 | |
---|
1155 | G4double piInMinusQ = std::pow(pi,-0.25) ; // 1.0/std::sqrt(std::sqrt(pi)) ?? |
---|
1156 | |
---|
1157 | G4int fNumber = (nHermite +1)/2 ; |
---|
1158 | G4double* fAbscissa = new G4double[fNumber] ; |
---|
1159 | G4double* fWeight = new G4double[fNumber] ; |
---|
1160 | |
---|
1161 | for(i=1;i<=fNumber;i++) |
---|
1162 | { |
---|
1163 | if(i == 1) |
---|
1164 | { |
---|
1165 | newton = std::sqrt((G4double)(2*nHermite + 1)) - |
---|
1166 | 1.85575001*std::pow((G4double)(2*nHermite + 1),-0.16666999) ; |
---|
1167 | } |
---|
1168 | else if(i == 2) |
---|
1169 | { |
---|
1170 | newton -= 1.14001*std::pow((G4double)nHermite,0.425999)/newton ; |
---|
1171 | } |
---|
1172 | else if(i == 3) |
---|
1173 | { |
---|
1174 | newton = 1.86002*newton - 0.86002*fAbscissa[0] ; |
---|
1175 | } |
---|
1176 | else if(i == 4) |
---|
1177 | { |
---|
1178 | newton = 1.91001*newton - 0.91001*fAbscissa[1] ; |
---|
1179 | } |
---|
1180 | else |
---|
1181 | { |
---|
1182 | newton = 2.0*newton - fAbscissa[i - 3] ; |
---|
1183 | } |
---|
1184 | for(k=1;k<=maxNumber;k++) |
---|
1185 | { |
---|
1186 | temp1 = piInMinusQ ; |
---|
1187 | temp2 = 0.0 ; |
---|
1188 | |
---|
1189 | for(j=1;j<=nHermite;j++) |
---|
1190 | { |
---|
1191 | temp3 = temp2 ; |
---|
1192 | temp2 = temp1 ; |
---|
1193 | temp1 = newton*std::sqrt(2.0/j)*temp2 - |
---|
1194 | std::sqrt(((G4double)(j - 1))/j)*temp3 ; |
---|
1195 | } |
---|
1196 | temp = std::sqrt((G4double)2*nHermite)*temp2 ; |
---|
1197 | newton1 = newton ; |
---|
1198 | newton = newton1 - temp1/temp ; |
---|
1199 | |
---|
1200 | if(std::fabs(newton - newton1) <= tolerance) |
---|
1201 | { |
---|
1202 | break ; |
---|
1203 | } |
---|
1204 | } |
---|
1205 | if(k > maxNumber) |
---|
1206 | { |
---|
1207 | G4Exception("G4Integrator<T,F>::Hermite(...)", "Error", |
---|
1208 | FatalException, "Too many (>12) iterations."); |
---|
1209 | } |
---|
1210 | fAbscissa[i-1] = newton ; |
---|
1211 | fWeight[i-1] = 2.0/(temp*temp) ; |
---|
1212 | } |
---|
1213 | |
---|
1214 | // |
---|
1215 | // Integral calculation |
---|
1216 | // |
---|
1217 | |
---|
1218 | for(i=0;i<fNumber;i++) |
---|
1219 | { |
---|
1220 | integral += fWeight[i]*( (*f)(fAbscissa[i]) + (*f)(-fAbscissa[i]) ) ; |
---|
1221 | } |
---|
1222 | delete[] fAbscissa; |
---|
1223 | delete[] fWeight; |
---|
1224 | return integral ; |
---|
1225 | } |
---|
1226 | |
---|
1227 | //////////////////////////////////////////////////////////////////////////// |
---|
1228 | // |
---|
1229 | // Method involving Jacobi polynomials |
---|
1230 | // |
---|
1231 | //////////////////////////////////////////////////////////////////////////// |
---|
1232 | // |
---|
1233 | // Gauss-Jacobi method for integration of ((1-x)^alpha)*((1+x)^beta)*f(x) |
---|
1234 | // from minus unit to plus unit . |
---|
1235 | // |
---|
1236 | |
---|
1237 | template <class T, class F> |
---|
1238 | G4double G4Integrator<T,F>::Jacobi( T& typeT, F f, G4double alpha, |
---|
1239 | G4double beta, G4int nJacobi) |
---|
1240 | { |
---|
1241 | const G4double tolerance = 1.0e-12 ; |
---|
1242 | const G4double maxNumber = 12 ; |
---|
1243 | G4int i, k, j ; |
---|
1244 | G4double alphaBeta, alphaReduced, betaReduced, root1=0., root2=0., root3=0. ; |
---|
1245 | G4double a, b, c, newton1, newton2, newton3, newton, temp, root=0., rootTemp ; |
---|
1246 | |
---|
1247 | G4int fNumber = nJacobi ; |
---|
1248 | G4double* fAbscissa = new G4double[fNumber] ; |
---|
1249 | G4double* fWeight = new G4double[fNumber] ; |
---|
1250 | |
---|
1251 | for (i=1;i<=nJacobi;i++) |
---|
1252 | { |
---|
1253 | if (i == 1) |
---|
1254 | { |
---|
1255 | alphaReduced = alpha/nJacobi ; |
---|
1256 | betaReduced = beta/nJacobi ; |
---|
1257 | root1 = (1.0+alpha)*(2.78002/(4.0+nJacobi*nJacobi)+ |
---|
1258 | 0.767999*alphaReduced/nJacobi) ; |
---|
1259 | root2 = 1.0+1.48*alphaReduced+0.96002*betaReduced + |
---|
1260 | 0.451998*alphaReduced*alphaReduced + |
---|
1261 | 0.83001*alphaReduced*betaReduced ; |
---|
1262 | root = 1.0-root1/root2 ; |
---|
1263 | } |
---|
1264 | else if (i == 2) |
---|
1265 | { |
---|
1266 | root1=(4.1002+alpha)/((1.0+alpha)*(1.0+0.155998*alpha)) ; |
---|
1267 | root2=1.0+0.06*(nJacobi-8.0)*(1.0+0.12*alpha)/nJacobi ; |
---|
1268 | root3=1.0+0.012002*beta*(1.0+0.24997*std::fabs(alpha))/nJacobi ; |
---|
1269 | root -= (1.0-root)*root1*root2*root3 ; |
---|
1270 | } |
---|
1271 | else if (i == 3) |
---|
1272 | { |
---|
1273 | root1=(1.67001+0.27998*alpha)/(1.0+0.37002*alpha) ; |
---|
1274 | root2=1.0+0.22*(nJacobi-8.0)/nJacobi ; |
---|
1275 | root3=1.0+8.0*beta/((6.28001+beta)*nJacobi*nJacobi) ; |
---|
1276 | root -= (fAbscissa[0]-root)*root1*root2*root3 ; |
---|
1277 | } |
---|
1278 | else if (i == nJacobi-1) |
---|
1279 | { |
---|
1280 | root1=(1.0+0.235002*beta)/(0.766001+0.118998*beta) ; |
---|
1281 | root2=1.0/(1.0+0.639002*(nJacobi-4.0)/(1.0+0.71001*(nJacobi-4.0))) ; |
---|
1282 | root3=1.0/(1.0+20.0*alpha/((7.5+alpha)*nJacobi*nJacobi)) ; |
---|
1283 | root += (root-fAbscissa[nJacobi-4])*root1*root2*root3 ; |
---|
1284 | } |
---|
1285 | else if (i == nJacobi) |
---|
1286 | { |
---|
1287 | root1 = (1.0+0.37002*beta)/(1.67001+0.27998*beta) ; |
---|
1288 | root2 = 1.0/(1.0+0.22*(nJacobi-8.0)/nJacobi) ; |
---|
1289 | root3 = 1.0/(1.0+8.0*alpha/((6.28002+alpha)*nJacobi*nJacobi)) ; |
---|
1290 | root += (root-fAbscissa[nJacobi-3])*root1*root2*root3 ; |
---|
1291 | } |
---|
1292 | else |
---|
1293 | { |
---|
1294 | root = 3.0*fAbscissa[i-2]-3.0*fAbscissa[i-3]+fAbscissa[i-4] ; |
---|
1295 | } |
---|
1296 | alphaBeta = alpha + beta ; |
---|
1297 | for (k=1;k<=maxNumber;k++) |
---|
1298 | { |
---|
1299 | temp = 2.0 + alphaBeta ; |
---|
1300 | newton1 = (alpha-beta+temp*root)/2.0 ; |
---|
1301 | newton2 = 1.0 ; |
---|
1302 | for (j=2;j<=nJacobi;j++) |
---|
1303 | { |
---|
1304 | newton3 = newton2 ; |
---|
1305 | newton2 = newton1 ; |
---|
1306 | temp = 2*j+alphaBeta ; |
---|
1307 | a = 2*j*(j+alphaBeta)*(temp-2.0) ; |
---|
1308 | b = (temp-1.0)*(alpha*alpha-beta*beta+temp*(temp-2.0)*root) ; |
---|
1309 | c = 2.0*(j-1+alpha)*(j-1+beta)*temp ; |
---|
1310 | newton1 = (b*newton2-c*newton3)/a ; |
---|
1311 | } |
---|
1312 | newton = (nJacobi*(alpha - beta - temp*root)*newton1 + |
---|
1313 | 2.0*(nJacobi + alpha)*(nJacobi + beta)*newton2)/ |
---|
1314 | (temp*(1.0 - root*root)) ; |
---|
1315 | rootTemp = root ; |
---|
1316 | root = rootTemp - newton1/newton ; |
---|
1317 | if (std::fabs(root-rootTemp) <= tolerance) |
---|
1318 | { |
---|
1319 | break ; |
---|
1320 | } |
---|
1321 | } |
---|
1322 | if (k > maxNumber) |
---|
1323 | { |
---|
1324 | G4Exception("G4Integrator<T,F>::Jacobi(T,F, ...)", "Error", |
---|
1325 | FatalException, "Too many (>12) iterations."); |
---|
1326 | } |
---|
1327 | fAbscissa[i-1] = root ; |
---|
1328 | fWeight[i-1] = std::exp(GammaLogarithm((G4double)(alpha+nJacobi)) + |
---|
1329 | GammaLogarithm((G4double)(beta+nJacobi)) - |
---|
1330 | GammaLogarithm((G4double)(nJacobi+1.0)) - |
---|
1331 | GammaLogarithm((G4double)(nJacobi + alphaBeta + 1.0))) |
---|
1332 | *temp*std::pow(2.0,alphaBeta)/(newton*newton2) ; |
---|
1333 | } |
---|
1334 | |
---|
1335 | // |
---|
1336 | // Calculation of the integral |
---|
1337 | // |
---|
1338 | |
---|
1339 | G4double integral = 0.0 ; |
---|
1340 | for(i=0;i<fNumber;i++) |
---|
1341 | { |
---|
1342 | integral += fWeight[i]*(typeT.*f)(fAbscissa[i]) ; |
---|
1343 | } |
---|
1344 | delete[] fAbscissa; |
---|
1345 | delete[] fWeight; |
---|
1346 | return integral ; |
---|
1347 | } |
---|
1348 | |
---|
1349 | |
---|
1350 | ///////////////////////////////////////////////////////////////////////// |
---|
1351 | // |
---|
1352 | // For use with 'this' pointer |
---|
1353 | |
---|
1354 | template <class T, class F> |
---|
1355 | G4double G4Integrator<T,F>::Jacobi( T* ptrT, F f, G4double alpha, |
---|
1356 | G4double beta, G4int n) |
---|
1357 | { |
---|
1358 | return Jacobi(*ptrT,f,alpha,beta,n) ; |
---|
1359 | } |
---|
1360 | |
---|
1361 | ///////////////////////////////////////////////////////////////////////// |
---|
1362 | // |
---|
1363 | // For use with global scope f |
---|
1364 | |
---|
1365 | template <class T, class F> |
---|
1366 | G4double G4Integrator<T,F>::Jacobi( G4double (*f)(G4double), G4double alpha, |
---|
1367 | G4double beta, G4int nJacobi) |
---|
1368 | { |
---|
1369 | const G4double tolerance = 1.0e-12 ; |
---|
1370 | const G4double maxNumber = 12 ; |
---|
1371 | G4int i, k, j ; |
---|
1372 | G4double alphaBeta, alphaReduced, betaReduced, root1=0., root2=0., root3=0. ; |
---|
1373 | G4double a, b, c, newton1, newton2, newton3, newton, temp, root=0., rootTemp ; |
---|
1374 | |
---|
1375 | G4int fNumber = nJacobi ; |
---|
1376 | G4double* fAbscissa = new G4double[fNumber] ; |
---|
1377 | G4double* fWeight = new G4double[fNumber] ; |
---|
1378 | |
---|
1379 | for (i=1;i<=nJacobi;i++) |
---|
1380 | { |
---|
1381 | if (i == 1) |
---|
1382 | { |
---|
1383 | alphaReduced = alpha/nJacobi ; |
---|
1384 | betaReduced = beta/nJacobi ; |
---|
1385 | root1 = (1.0+alpha)*(2.78002/(4.0+nJacobi*nJacobi)+ |
---|
1386 | 0.767999*alphaReduced/nJacobi) ; |
---|
1387 | root2 = 1.0+1.48*alphaReduced+0.96002*betaReduced + |
---|
1388 | 0.451998*alphaReduced*alphaReduced + |
---|
1389 | 0.83001*alphaReduced*betaReduced ; |
---|
1390 | root = 1.0-root1/root2 ; |
---|
1391 | } |
---|
1392 | else if (i == 2) |
---|
1393 | { |
---|
1394 | root1=(4.1002+alpha)/((1.0+alpha)*(1.0+0.155998*alpha)) ; |
---|
1395 | root2=1.0+0.06*(nJacobi-8.0)*(1.0+0.12*alpha)/nJacobi ; |
---|
1396 | root3=1.0+0.012002*beta*(1.0+0.24997*std::fabs(alpha))/nJacobi ; |
---|
1397 | root -= (1.0-root)*root1*root2*root3 ; |
---|
1398 | } |
---|
1399 | else if (i == 3) |
---|
1400 | { |
---|
1401 | root1=(1.67001+0.27998*alpha)/(1.0+0.37002*alpha) ; |
---|
1402 | root2=1.0+0.22*(nJacobi-8.0)/nJacobi ; |
---|
1403 | root3=1.0+8.0*beta/((6.28001+beta)*nJacobi*nJacobi) ; |
---|
1404 | root -= (fAbscissa[0]-root)*root1*root2*root3 ; |
---|
1405 | } |
---|
1406 | else if (i == nJacobi-1) |
---|
1407 | { |
---|
1408 | root1=(1.0+0.235002*beta)/(0.766001+0.118998*beta) ; |
---|
1409 | root2=1.0/(1.0+0.639002*(nJacobi-4.0)/(1.0+0.71001*(nJacobi-4.0))) ; |
---|
1410 | root3=1.0/(1.0+20.0*alpha/((7.5+alpha)*nJacobi*nJacobi)) ; |
---|
1411 | root += (root-fAbscissa[nJacobi-4])*root1*root2*root3 ; |
---|
1412 | } |
---|
1413 | else if (i == nJacobi) |
---|
1414 | { |
---|
1415 | root1 = (1.0+0.37002*beta)/(1.67001+0.27998*beta) ; |
---|
1416 | root2 = 1.0/(1.0+0.22*(nJacobi-8.0)/nJacobi) ; |
---|
1417 | root3 = 1.0/(1.0+8.0*alpha/((6.28002+alpha)*nJacobi*nJacobi)) ; |
---|
1418 | root += (root-fAbscissa[nJacobi-3])*root1*root2*root3 ; |
---|
1419 | } |
---|
1420 | else |
---|
1421 | { |
---|
1422 | root = 3.0*fAbscissa[i-2]-3.0*fAbscissa[i-3]+fAbscissa[i-4] ; |
---|
1423 | } |
---|
1424 | alphaBeta = alpha + beta ; |
---|
1425 | for (k=1;k<=maxNumber;k++) |
---|
1426 | { |
---|
1427 | temp = 2.0 + alphaBeta ; |
---|
1428 | newton1 = (alpha-beta+temp*root)/2.0 ; |
---|
1429 | newton2 = 1.0 ; |
---|
1430 | for (j=2;j<=nJacobi;j++) |
---|
1431 | { |
---|
1432 | newton3 = newton2 ; |
---|
1433 | newton2 = newton1 ; |
---|
1434 | temp = 2*j+alphaBeta ; |
---|
1435 | a = 2*j*(j+alphaBeta)*(temp-2.0) ; |
---|
1436 | b = (temp-1.0)*(alpha*alpha-beta*beta+temp*(temp-2.0)*root) ; |
---|
1437 | c = 2.0*(j-1+alpha)*(j-1+beta)*temp ; |
---|
1438 | newton1 = (b*newton2-c*newton3)/a ; |
---|
1439 | } |
---|
1440 | newton = (nJacobi*(alpha - beta - temp*root)*newton1 + |
---|
1441 | 2.0*(nJacobi + alpha)*(nJacobi + beta)*newton2) / |
---|
1442 | (temp*(1.0 - root*root)) ; |
---|
1443 | rootTemp = root ; |
---|
1444 | root = rootTemp - newton1/newton ; |
---|
1445 | if (std::fabs(root-rootTemp) <= tolerance) |
---|
1446 | { |
---|
1447 | break ; |
---|
1448 | } |
---|
1449 | } |
---|
1450 | if (k > maxNumber) |
---|
1451 | { |
---|
1452 | G4Exception("G4Integrator<T,F>::Jacobi(...)", "Error", |
---|
1453 | FatalException, "Too many (>12) iterations."); |
---|
1454 | } |
---|
1455 | fAbscissa[i-1] = root ; |
---|
1456 | fWeight[i-1] = |
---|
1457 | std::exp(GammaLogarithm((G4double)(alpha+nJacobi)) + |
---|
1458 | GammaLogarithm((G4double)(beta+nJacobi)) - |
---|
1459 | GammaLogarithm((G4double)(nJacobi+1.0)) - |
---|
1460 | GammaLogarithm((G4double)(nJacobi + alphaBeta + 1.0))) |
---|
1461 | *temp*std::pow(2.0,alphaBeta)/(newton*newton2); |
---|
1462 | } |
---|
1463 | |
---|
1464 | // |
---|
1465 | // Calculation of the integral |
---|
1466 | // |
---|
1467 | |
---|
1468 | G4double integral = 0.0 ; |
---|
1469 | for(i=0;i<fNumber;i++) |
---|
1470 | { |
---|
1471 | integral += fWeight[i]*(*f)(fAbscissa[i]) ; |
---|
1472 | } |
---|
1473 | delete[] fAbscissa; |
---|
1474 | delete[] fWeight; |
---|
1475 | return integral ; |
---|
1476 | } |
---|
1477 | |
---|
1478 | // |
---|
1479 | // |
---|
1480 | /////////////////////////////////////////////////////////////////// |
---|