source: trunk/source/parameterisations/gflash/src/Gamma.cc@ 1057

Last change on this file since 1057 was 850, checked in by garnier, 17 years ago

geant4.8.2 beta

File size: 5.8 KB
Line 
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// $Id: Gamma.cc,v 1.6 2006/06/29 19:14:28 gunter Exp $
27// GEANT4 tag $Name: HEAD $
28//
29//
30// ------------------------------------------------------------
31// GEANT 4 class implementation
32// ------------------------------------------------------------
33
34#include <cmath>
35#include <string.h>
36#include "Gamma.hh"
37
38MyGamma::MyGamma(){}
39
40MyGamma::~MyGamma(){}
41
42//____________________________________________________________________________
43double MyGamma::Gamma(double z)
44{
45 // Computation of gamma(z) for all z>0.
46 //
47 // The algorithm is based on the article by C.Lanczos [1] as denoted in
48 // Numerical Recipes 2nd ed. on p. 207 (W.H.Press et al.).
49 //
50 // [1] C.Lanczos, SIAM Journal of Numerical Analysis B1 (1964), 86.
51 //
52 //--- Nve 14-nov-1998 UU-SAP Utrecht
53
54 if (z<=0) return 0;
55
56 double v = LnGamma(z);
57 return std::exp(v);
58}
59
60//____________________________________________________________________________
61double MyGamma::Gamma(double a,double x)
62{
63 // Computation of the incomplete gamma function P(a,x)
64 //
65 // The algorithm is based on the formulas and code as denoted in
66 // Numerical Recipes 2nd ed. on p. 210-212 (W.H.Press et al.).
67 //
68 //--- Nve 14-nov-1998 UU-SAP Utrecht
69
70 if (a <= 0 || x <= 0) return 0;
71
72 if (x < (a+1)) return GamSer(a,x);
73 else return GamCf(a,x);
74}
75
76//____________________________________________________________________________
77double MyGamma::GamCf(double a,double x)
78{
79 // Computation of the incomplete gamma function P(a,x)
80 // via its continued fraction representation.
81 //
82 // The algorithm is based on the formulas and code as denoted in
83 // Numerical Recipes 2nd ed. on p. 210-212 (W.H.Press et al.).
84 //
85 //--- Nve 14-nov-1998 UU-SAP Utrecht
86
87 int itmax = 100; // Maximum number of iterations
88 double eps = 3.e-7; // Relative accuracy
89 double fpmin = 1.e-30; // Smallest double value allowed here
90
91 if (a <= 0 || x <= 0) return 0;
92
93 double gln = LnGamma(a);
94 double b = x+1-a;
95 double c = 1/fpmin;
96 double d = 1/b;
97 double h = d;
98 double an,del;
99 for (int i=1; i<=itmax; i++) {
100 an = double(-i)*(double(i)-a);
101 b += 2;
102 d = an*d+b;
103 if (Abs(d) < fpmin) d = fpmin;
104 c = b+an/c;
105 if (Abs(c) < fpmin) c = fpmin;
106 d = 1/d;
107 del = d*c;
108 h = h*del;
109 if (Abs(del-1) < eps) break;
110 //if (i==itmax) cout << "*GamCf(a,x)* a too large or itmax too small" << endl;
111 }
112 double v = Exp(-x+a*Log(x)-gln)*h;
113 return (1-v);
114}
115
116//____________________________________________________________________________
117double MyGamma::GamSer(double a,double x)
118{
119 // Computation of the incomplete gamma function P(a,x)
120 // via its series representation.
121 //
122 // The algorithm is based on the formulas and code as denoted in
123 // Numerical Recipes 2nd ed. on p. 210-212 (W.H.Press et al.).
124 //
125 //--- Nve 14-nov-1998 UU-SAP Utrecht
126
127 int itmax = 100; // Maximum number of iterations
128 double eps = 3.e-7; // Relative accuracy
129
130 if (a <= 0 || x <= 0) return 0;
131
132 double gln = LnGamma(a);
133 double ap = a;
134 double sum = 1/a;
135 double del = sum;
136 for (int n=1; n<=itmax; n++) {
137 ap += 1;
138 del = del*x/ap;
139 sum += del;
140 if (MyGamma::Abs(del) < Abs(sum*eps)) break;
141 //if (n==itmax) cout << "*GamSer(a,x)* a too large or itmax too small" << endl;
142 }
143 double v = sum*Exp(-x+a*Log(x)-gln);
144 return v;
145}
146
147
148double MyGamma::LnGamma(double z)
149{
150 // Computation of ln[gamma(z)] for all z>0.
151 //
152 // The algorithm is based on the article by C.Lanczos [1] as denoted in
153 // Numerical Recipes 2nd ed. on p. 207 (W.H.Press et al.).
154 //
155 // [1] C.Lanczos, SIAM Journal of Numerical Analysis B1 (1964), 86.
156 //
157 // The accuracy of the result is better than 2e-10.
158 //
159 //--- Nve 14-nov-1998 UU-SAP Utrecht
160
161 if (z<=0) return 0;
162
163 // Coefficients for the series expansion
164 double c[7] = { 2.5066282746310005, 76.18009172947146, -86.50532032941677
165 ,24.01409824083091, -1.231739572450155, 0.1208650973866179e-2
166 ,-0.5395239384953e-5};
167
168 double x = z;
169 double y = x;
170 double tmp = x+5.5;
171 tmp = (x+0.5)*Log(tmp)-tmp;
172 double ser = 1.000000000190015;
173 for (int i=1; i<7; i++) {
174 y += 1;
175 ser += c[i]/y;
176 }
177 double v = tmp+Log(c[0]*ser/x);
178 return v;
179}
Note: See TracBrowser for help on using the repository browser.