Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
nf_specialFunctions.h File Reference
#include <math.h>
#include <float.h>
#include <nf_utilities.h>
Include dependency graph for nf_specialFunctions.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

double nf_polevl (double x, double coef[], int N)
 
double nf_p1evl (double x, double coef[], int N)
 
double nf_exponentialIntegral (int n, double x, nfu_status *status)
 
double nf_gammaFunction (double x, nfu_status *status)
 
double nf_logGammaFunction (double x, nfu_status *status)
 
double nf_incompleteGammaFunction (double a, double x, nfu_status *status)
 
double nf_incompleteGammaFunctionComplementary (double a, double x, nfu_status *status)
 
double nf_amc_log_factorial (int)
 
double nf_amc_factorial (int)
 
double nf_amc_wigner_3j (int, int, int, int, int, int)
 
double nf_amc_wigner_6j (int, int, int, int, int, int)
 
double nf_amc_wigner_9j (int, int, int, int, int, int, int, int, int)
 
double nf_amc_racah (int, int, int, int, int, int)
 
double nf_amc_clebsh_gordan (int, int, int, int, int)
 
double nf_amc_z_coefficient (int, int, int, int, int, int)
 
double nf_amc_zbar_coefficient (int, int, int, int, int, int)
 
double nf_amc_reduced_matrix_element (int, int, int, int, int, int, int)
 

Function Documentation

double nf_amc_clebsh_gordan ( int  ,
int  ,
int  ,
int  ,
int   
)

Definition at line 286 of file nf_angularMomentumCoupling.cc.

286  {
287 /*
288 * Clebsh-Gordan coefficient
289 * = <j1,j2,m1,m2|j3,m1+m2>
290 * = (-)^(j1-j2+m1+m2) * std::sqrt(2*j3+1) * / j1 j2 j3 \
291 * \ m1 m2 -m1-m2 /
292 *
293 * Note: Last value m3 is preset to m1+m2. Any other value will evaluate to 0.0.
294 */
295 
296  int m3, x1, x2, x3, y1, y2, y3;
297  double cg = 0.0;
298 
299  if ( j1 < 0 || j2 < 0 || j3 < 0) return( 0.0 );
300  if ( j1 + j2 + j3 > 2 * MAX_FACTORIAL ) return( INFINITY );
301 
302  m3 = m1 + m2;
303 
304  if ( ( x1 = ( j1 + m1 ) / 2 + 1 ) <= 0 ) return( 0.0 );
305  if ( ( x2 = ( j2 + m2 ) / 2 + 1 ) <= 0 ) return( 0.0 );
306  if ( ( x3 = ( j3 - m3 ) / 2 + 1 ) <= 0 ) return( 0.0 );
307 
308  if ( ( y1 = x1 - m1 ) <= 0 ) return( 0.0 );
309  if ( ( y2 = x2 - m2 ) <= 0 ) return( 0.0 );
310  if ( ( y3 = x3 + m3 ) <= 0 ) return( 0.0 );
311 
312  if ( j3 == 0 ){
313  if ( j1 == j2 ) cg = ( 1.0 / std::sqrt( (double)j1 + 1.0 ) * ( ( y1 % 2 == 0 ) ? -1:1 ) );
314  }
315  else if ( (j1 == 0 || j2 == 0 ) ){
316  if ( ( j1 + j2 ) == j3 ) cg = 1.0;
317  }
318  else {
319  if( m3 == 0 && std::abs( m1 ) <= 1 ){
320  if( m1 == 0 ) cg = cg1( x1, x2, x3 );
321  else cg = cg2( x1 + y1 - y2, x3 - 1, x1 + x2 - 2, x1 - y2, j1, j2, j3, m2 );
322  }
323  else if ( m2 == 0 && std::abs( m1 ) <=1 ){
324  cg = cg2( x1 - y2 + y3, x2 - 1, x1 + x3 - 2, x3 - y1, j1, j3, j3, m1 );
325  }
326  else if ( m1 == 0 && std::abs( m3 ) <= 1 ){
327  cg = cg2( x1, x1 - 1, x2 + x3 - 2, x2 - y3, j2, j3, j3, -m3 );
328  }
329  else cg = cg3( x1, x2, x3, y1, y2, y3 );
330  }
331 
332  return( cg );
333 }
static double cg3(int, int, int, int, int, int)
static constexpr double m3
Definition: G4SIunits.hh:131
static double cg1(int, int, int)
static const int MAX_FACTORIAL
static double cg2(int, int, int, int, int, int, int, int)
static constexpr double m2
Definition: G4SIunits.hh:130

Here is the call graph for this function:

Here is the caller graph for this function:

double nf_amc_factorial ( int  )

Definition at line 96 of file nf_angularMomentumCoupling.cc.

96  {
97 /*
98 * returns n! for pre-computed table. INFINITY is return if n is negative or too large.
99 */
100  return G4Exp( nf_amc_log_factorial( n ) );
101 }
double nf_amc_log_factorial(int)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183

Here is the call graph for this function:

double nf_amc_log_factorial ( int  )

Definition at line 85 of file nf_angularMomentumCoupling.cc.

85  {
86 /*
87 * returns ln( n! ).
88 */
89  if( n > MAX_FACTORIAL ) return( INFINITY );
90  if( n < 0 ) return( INFINITY );
91  return nf_amc_log_fact[n];
92 }
static const int MAX_FACTORIAL
static const double nf_amc_log_fact[]

Here is the caller graph for this function:

double nf_amc_racah ( int  ,
int  ,
int  ,
int  ,
int  ,
int   
)

Definition at line 251 of file nf_angularMomentumCoupling.cc.

251  {
252 /*
253 * Racah coefficient definition in Edmonds (AR Edmonds, "Angular Momentum in Quantum Mechanics", Princeton (1980) is
254 * W(j1, j2, l2, l1 ; j3, l3) = (-1)^(j1+j2+l1+l2) * { j1 j2 j3 }
255 * { l1 l2 l3 }
256 * The call signature of W(...) appears jumbled, but hey, that's the convention.
257 *
258 * This convention is exactly that used by Blatt-Biedenharn (Rev. Mod. Phys. 24, 258 (1952)) too
259 */
260 
261  double sig;
262 
263  sig = ( ( ( j1 + j2 + l1 + l2 ) % 4 == 0 ) ? 1.0 : -1.0 );
264  return sig * nf_amc_wigner_6j( j1, j2, j3, l1, l2, l3 );
265 }
double nf_amc_wigner_6j(int, int, int, int, int, int)

Here is the call graph for this function:

Here is the caller graph for this function:

double nf_amc_reduced_matrix_element ( int  ,
int  ,
int  ,
int  ,
int  ,
int  ,
int   
)

Definition at line 471 of file nf_angularMomentumCoupling.cc.

471  {
472 /*
473 * Reduced Matrix Element for Tensor Operator
474 * = < l1j1 || T(YL,sigma_S)J || l0j0 >
475 *
476 * M.B.Johnson, L.W.Owen, G.R.Satchler
477 * Phys. Rev. 142, 748 (1966)
478 * Note: definition differs from JOS by the factor sqrt(2j1+1)
479 */
480  int llt;
481  double x1, x2, x3, reduced_mat, clebsh_gordan;
482 
483  if ( parity( lt ) != parity( l0 ) * parity( l1 ) ) return( 0.0 );
484  if ( std::abs( l0 - l1 ) > lt || ( l0 + l1 ) < lt ) return( 0.0 );
485  if ( std::abs( ( j0 - j1 ) / 2 ) > jt || ( ( j0 + j1 ) / 2 ) < jt ) return( 0.0 );
486 
487  llt = 2 * lt;
488  jt *= 2;
489  st *= 2;
490 
491  if( ( clebsh_gordan = nf_amc_clebsh_gordan( j1, j0, 1, -1, jt ) ) == INFINITY ) return( INFINITY );
492 
493  reduced_mat = 1.0 / std::sqrt( 4 * M_PI ) * clebsh_gordan / std::sqrt( jt + 1.0 ) /* BRB jt + 1 <= 0? */
494  * std::sqrt( ( j0 + 1.0 ) * ( j1 + 1.0 ) * ( llt + 1.0 ) )
495  * parity( ( j1 - j0 ) / 2 ) * parity( ( -l0 + l1 + lt ) / 2 ) * parity( ( j0 - 1 ) / 2 );
496 
497  if( st == 2 ){
498  x1 = ( l0 - j0 / 2.0 ) * ( j0 + 1.0 );
499  x2 = ( l1 - j1 / 2.0 ) * ( j1 + 1.0 );
500  if ( jt == llt ){
501  x3 = ( lt == 0 ) ? 0 : ( x1 - x2 ) / std::sqrt( lt * ( lt + 1.0 ) );
502  }
503  else if ( jt == ( llt - st ) ){
504  x3 = ( lt == 0 ) ? 0 : -( lt + x1 + x2 ) / std::sqrt( lt * ( 2.0 * lt + 1.0 ) );
505  }
506  else if ( jt == ( llt + st ) ){
507  x3 = ( lt + 1 - x1 - x2 ) / std::sqrt( ( 2.0 * lt + 1.0 ) * ( lt + 1.0 ) );
508  }
509  else{
510  x3 = 1.0;
511  }
512  }
513  else x3 = 1.0;
514  reduced_mat *= x3;
515 
516  return( reduced_mat );
517 }
#define M_PI
Definition: SbMath.h:34
static int parity(int x)
double nf_amc_clebsh_gordan(int, int, int, int, int)

Here is the call graph for this function:

double nf_amc_wigner_3j ( int  ,
int  ,
int  ,
int  ,
int  ,
int   
)

Definition at line 105 of file nf_angularMomentumCoupling.cc.

105  {
106 /*
107 * Wigner's 3J symbol (similar to Clebsh-Gordan)
108 * = / j1 j2 j3 \
109 * \ j4 j5 j6 /
110 */
111  double cg;
112 
113  if( ( j4 + j5 + j6 ) != 0 ) return( 0.0 );
114  if( ( cg = nf_amc_clebsh_gordan( j1, j2, j4, j5, j3 ) ) == 0.0 ) return ( 0.0 );
115  if( cg == INFINITY ) return( cg );
116  return( ( ( ( j1 - j2 - j6 ) % 4 == 0 ) ? 1.0 : -1.0 ) * cg / std::sqrt( j3 + 1.0 ) ); /* BRB j3 + 1 <= 0? */
117 }
double nf_amc_clebsh_gordan(int, int, int, int, int)

Here is the call graph for this function:

double nf_amc_wigner_6j ( int  ,
int  ,
int  ,
int  ,
int  ,
int   
)

Definition at line 121 of file nf_angularMomentumCoupling.cc.

121  {
122 /*
123 * Wigner's 6J symbol (similar to Racah)
124 * = { j1 j2 j3 }
125 * { j4 j5 j6 }
126 */
127  int i, x[6];
128 
129  x[0] = j1; x[1] = j2; x[2] = j3; x[3] = j4; x[4] = j5; x[5] = j6;
130  for( i = 0; i < 6; i++ ) if ( x[i] == 0 ) return( w6j0( i, x ) );
131 
132  return( w6j1( x ) );
133 }
static double w6j0(int, int *)
static double w6j1(int *)

Here is the call graph for this function:

Here is the caller graph for this function:

double nf_amc_wigner_9j ( int  ,
int  ,
int  ,
int  ,
int  ,
int  ,
int  ,
int  ,
int   
)

Definition at line 224 of file nf_angularMomentumCoupling.cc.

224  {
225 /*
226 * Wigner's 9J symbol
227 * / j1 j2 j3 \
228 * = | j4 j5 j6 |
229 * \ j7 j8 j9 /
230 *
231 */
232  int i, i0, i1;
233  double rac;
234 
235  i0 = max3( std::abs( j1 - j9 ), std::abs( j2 - j6 ), std::abs( j4 - j8 ) );
236  i1 = min3( ( j1 + j9 ), ( j2 + j6 ), ( j4 + j8 ) );
237 
238  rac = 0.0;
239  for ( i = i0; i <= i1; i += 2 ){
240  rac += nf_amc_racah( j1, j4, j9, j8, j7, i )
241  * nf_amc_racah( j2, j5, i, j4, j8, j6 )
242  * nf_amc_racah( j9, i, j3, j2, j1, j6 ) * ( i + 1 );
243  if( rac == INFINITY ) return( INFINITY );
244  }
245 
246  return( ( ( (int)( ( j1 + j3 + j5 + j8 ) / 2 + j2 + j4 + j9 ) % 4 == 0 ) ? 1.0 : -1.0 ) * rac );
247 }
static int min3(int a, int b, int c)
double nf_amc_racah(int, int, int, int, int, int)
static int max3(int a, int b, int c)

Here is the call graph for this function:

double nf_amc_z_coefficient ( int  ,
int  ,
int  ,
int  ,
int  ,
int   
)

Definition at line 435 of file nf_angularMomentumCoupling.cc.

435  {
436 /*
437 * Biedenharn's Z-coefficient coefficient
438 * = Z(l1 j1 l2 j2 | S L )
439 */
440  double z, clebsh_gordan = nf_amc_clebsh_gordan( l1, l2, 0, 0, ll ), racah = nf_amc_racah( l1, j1, l2, j2, s, ll );
441 
442  if( ( clebsh_gordan == INFINITY ) || ( racah == INFINITY ) ) return( INFINITY );
443  z = ( ( ( -l1 + l2 + ll ) % 8 == 0 ) ? 1.0 : -1.0 )
444  * std::sqrt( l1 + 1.0 ) * std::sqrt( l2 + 1.0 ) * std::sqrt( j1 + 1.0 ) * std::sqrt( j2 + 1.0 ) * clebsh_gordan * racah;
445 
446  return( z );
447 }
const XML_Char * s
Definition: expat.h:262
double nf_amc_racah(int, int, int, int, int, int)
double nf_amc_clebsh_gordan(int, int, int, int, int)

Here is the call graph for this function:

double nf_amc_zbar_coefficient ( int  ,
int  ,
int  ,
int  ,
int  ,
int   
)

Definition at line 451 of file nf_angularMomentumCoupling.cc.

451  {
452 /*
453 * Lane & Thomas's Zbar-coefficient coefficient
454 * = Zbar(l1 j1 l2 j2 | S L )
455 * = (-i)^( -l1 + l2 + ll ) * Z(l1 j1 l2 j2 | S L )
456 *
457 * Lane & Thomas Rev. Mod. Phys. 30, 257-353 (1958).
458 * Note, Lane & Thomas define this because they did not like the different phase convention in Blatt & Biedenharn's Z coefficient. They changed it to get better time-reversal behavior.
459 * Froehner uses Lane & Thomas convention as does T. Kawano.
460 */
461  double zbar, clebsh_gordan = nf_amc_clebsh_gordan( l1, l2, 0, 0, ll ), racah = nf_amc_racah( l1, j1, l2, j2, s, ll );
462 
463  if( ( clebsh_gordan == INFINITY ) || ( racah == INFINITY ) ) return( INFINITY );
464  zbar = std::sqrt( l1 + 1.0 ) * std::sqrt( l2 + 1.0 ) * std::sqrt( j1 + 1.0 ) * std::sqrt( j2 + 1.0 ) * clebsh_gordan * racah;
465 
466  return( zbar );
467 }
const XML_Char * s
Definition: expat.h:262
double nf_amc_racah(int, int, int, int, int, int)
double nf_amc_clebsh_gordan(int, int, int, int, int)

Here is the call graph for this function:

double nf_exponentialIntegral ( int  n,
double  x,
nfu_status status 
)

Definition at line 28 of file nf_exponentialIntegral.cc.

28  {
29 
30  int i, ii, nm1;
31  double a, b, c, d, del, fact, h, psi;
32  double ans = 0.0;
33 
34  *status = nfu_badInput;
35  if( !isfinite( x ) ) return( x );
36  *status = nfu_Okay;
37 
38  nm1 = n - 1;
39  if( ( n < 0 ) || ( x < 0.0 ) || ( ( x == 0.0 ) && ( ( n == 0 ) || ( n == 1 ) ) ) ) {
40  *status = nfu_badInput; }
41  else {
42  if( n == 0 ) {
43  ans = G4Exp( -x ) / x; } /* Special case */
44  else if( x == 0.0 ) {
45  ans = 1.0 / nm1; } /* Another special case */
46  else if( x > 1.0 ) { /* Lentz's algorithm */
47  b = x + n;
48  c = 1.0 / FPMIN;
49  d = 1.0 / b;
50  h = d;
51  for( i = 1; i <= MAXIT; i++ ) {
52  a = -i * ( nm1 + i );
53  b += 2.0;
54  d = 1.0 / ( a * d + b ); /* Denominators cannot be zero */
55  c = b + a / c;
56  del = c * d;
57  h *= del;
58  if( fabs( del - 1.0 ) < EPS ) return( h * G4Exp( -x ) );
59  }
60  *status = nfu_failedToConverge; }
61  else {
62  ans = ( nm1 != 0 ) ? 1.0 / nm1 : -G4Log(x) - EULER; /* Set first term */
63  fact = 1.0;
64  for( i = 1; i <= MAXIT; i++ ) {
65  fact *= -x / i;
66  if( i != nm1 ) {
67  del = -fact / ( i - nm1 ); }
68  else {
69  psi = -EULER; /* Compute psi(n) */
70  for( ii = 1; ii <= nm1; ii++ ) psi += 1.0 / ii;
71  del = fact * ( -G4Log( x ) + psi );
72  }
73  ans += del;
74  if( fabs( del ) < fabs( ans ) * EPS ) return( ans );
75  }
76  *status = nfu_failedToConverge;
77  }
78  }
79  return( ans );
80 }
#define FPMIN
#define MAXIT
#define isfinite
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
#define EULER
#define EPS

Here is the call graph for this function:

Here is the caller graph for this function:

double nf_gammaFunction ( double  x,
nfu_status status 
)

Definition at line 126 of file nf_gammaFunctions.cc.

126  {
127 
128  double p, q, z;
129  int i, sgngam = 1;
130 
131  *status = nfu_badInput;
132  if( !isfinite( x ) ) return( x );
133  *status = nfu_Okay;
134 
135  q = fabs( x );
136 
137  if( q > 33.0 ) {
138  if( x < 0.0 ) {
139  p = floor( q );
140  if( p == q ) goto goverf;
141  i = (int) p;
142  if( ( i & 1 ) == 0 ) sgngam = -1;
143  z = q - p;
144  if( z > 0.5 ) {
145  p += 1.0;
146  z = q - p;
147  }
148  z = q * sin( M_PI * z );
149  if( z == 0.0 ) goto goverf;
150  z = M_PI / ( fabs( z ) * stirf( q, status ) );
151  }
152  else {
153  z = stirf( x, status );
154  }
155  return( sgngam * z );
156  }
157 
158  z = 1.0;
159  while( x >= 3.0 ) {
160  x -= 1.0;
161  z *= x;
162  } // Loop checking, 11.06.2015, T. Koi
163 
164  while( x < 0.0 ) {
165  if( x > -1.E-9 ) goto small;
166  z /= x;
167  x += 1.0;
168  } // Loop checking, 11.06.2015, T. Koi
169 
170  while( x < 2.0 ) {
171  if( x < 1.e-9 ) goto small;
172  z /= x;
173  x += 1.0;
174  } // Loop checking, 11.06.2015, T. Koi
175 
176  if( x == 2.0 ) return( z );
177 
178  x -= 2.0;
179  p = nf_polevl( x, P, 6 );
180  q = nf_polevl( x, Q, 7 );
181  return( z * p / q );
182 
183 small:
184  if( x == 0.0 ) goto goverf;
185  return( z / ( ( 1.0 + 0.5772156649015329 * x ) * x ) );
186 
187 goverf:
188  return( sgngam * DBL_MAX );
189 }
#define M_PI
Definition: SbMath.h:34
const char * p
Definition: xmltok.h:285
double nf_polevl(double x, double coef[], int N)
Definition: nf_polevl.cc:46
static double Q[]
static double stirf(double x, nfu_status *status)
static double P[]
typedef int(XMLCALL *XML_NotStandaloneHandler)(void *userData)
#define isfinite
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:

Here is the caller graph for this function:

double nf_incompleteGammaFunction ( double  a,
double  x,
nfu_status status 
)

Definition at line 155 of file nf_incompleteGammaFunctions.cc.

155  {
156 /* left tail of incomplete gamma function:
157 *
158 * inf. k
159 * a -x - x
160 * x e > ----------
161 * - -
162 * k=0 | (a+k+1)
163 */
164  double ans, ax, c, r;
165 
166  *status = nfu_badInput;
167  if( !isfinite( x ) ) return( x );
168  *status = nfu_Okay;
169 
170  if( ( x <= 0 ) || ( a <= 0 ) ) return( 0.0 );
171  if( ( x > 1.0 ) && ( x > a ) ) return( nf_gammaFunction( a, status ) - nf_incompleteGammaFunctionComplementary( a, x, status ) );
172 
173  ax = G4Exp( a * G4Log( x ) - x ); /* Compute x**a * exp(-x) */
174  if( ax == 0. ) return( 0.0 );
175 
176  r = a; /* power series */
177  c = 1.0;
178  ans = 1.0;
179  do {
180  r += 1.0;
181  c *= x / r;
182  ans += c;
183  } while( c > ans * DBL_EPSILON ); // Loop checking, 11.06.2015, T. Koi
184 
185  return( ans * ax / a );
186 }
#define DBL_EPSILON
Definition: templates.hh:87
#define isfinite
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
double nf_incompleteGammaFunctionComplementary(double a, double x, nfu_status *status)
double nf_gammaFunction(double x, nfu_status *status)

Here is the call graph for this function:

Here is the caller graph for this function:

double nf_incompleteGammaFunctionComplementary ( double  a,
double  x,
nfu_status status 
)

Definition at line 88 of file nf_incompleteGammaFunctions.cc.

88  {
89 
90  double ans, ax, c, yc, r, t, y, z;
91  double pk, pkm1, pkm2, qk, qkm1, qkm2;
92 
93  *status = nfu_badInput;
94  if( !isfinite( x ) ) return( x );
95  *status = nfu_Okay;
96 
97  if( ( x <= 0 ) || ( a <= 0 ) ) return( 1.0 );
98  if( ( x < 1.0 ) || ( x < a ) ) return( nf_gammaFunction( a, status ) - nf_incompleteGammaFunction( a, x, status ) );
99 
100  ax = G4Exp( a * G4Log( x ) - x );
101  if( ax == 0. ) return( 0.0 );
102 
103  if( x < 10000. ) {
104  y = 1.0 - a; /* continued fraction */
105  z = x + y + 1.0;
106  c = 0.0;
107  pkm2 = 1.0;
108  qkm2 = x;
109  pkm1 = x + 1.0;
110  qkm1 = z * x;
111  ans = pkm1 / qkm1;
112 
113  do {
114  c += 1.0;
115  y += 1.0;
116  z += 2.0;
117  yc = y * c;
118  pk = pkm1 * z - pkm2 * yc;
119  qk = qkm1 * z - qkm2 * yc;
120  if( qk != 0 ) {
121  r = pk / qk;
122  t = fabs( ( ans - r ) / r );
123  ans = r; }
124  else {
125  t = 1.0;
126  }
127  pkm2 = pkm1;
128  pkm1 = pk;
129  qkm2 = qkm1;
130  qkm1 = qk;
131  if( fabs( pk ) > big ) {
132  pkm2 *= biginv;
133  pkm1 *= biginv;
134  qkm2 *= biginv;
135  qkm1 *= biginv;
136  }
137  } while( t > DBL_EPSILON ); } // Loop checking, 11.06.2015, T. Koi
138  else { /* Asymptotic expansion. */
139  y = 1. / x;
140  r = a;
141  c = 1.;
142  ans = 1.;
143  do {
144  a -= 1.;
145  c *= a * y;
146  ans += c;
147  } while( fabs( c ) > 100 * ans * DBL_EPSILON ); // Loop checking, 11.06.2015, T. Koi
148  }
149 
150  return( ans * ax );
151 }
static double biginv
#define DBL_EPSILON
Definition: templates.hh:87
#define isfinite
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
double nf_incompleteGammaFunction(double a, double x, nfu_status *status)
static double big
double nf_gammaFunction(double x, nfu_status *status)

Here is the call graph for this function:

Here is the caller graph for this function:

double nf_logGammaFunction ( double  x,
nfu_status status 
)

Definition at line 206 of file nf_gammaFunctions.cc.

206  {
207 /* Logarithm of gamma function */
208 
209  int sgngam;
210 
211  *status = nfu_badInput;
212  if( !isfinite( x ) ) return( x );
213  *status = nfu_Okay;
214  return( lgam( x, &sgngam, status ) );
215 }
#define isfinite
static double lgam(double x, int *sgngam, nfu_status *status)

Here is the call graph for this function:

double nf_p1evl ( double  x,
double  coef[],
int  N 
)

Definition at line 67 of file nf_polevl.cc.

67  {
68 
69  double ans;
70  double *p;
71  int i;
72 
73  p = coef;
74  ans = x + *p++;
75  i = N-1;
76 
77  do {
78  ans = ans * x + *p++; }
79  while( --i ); // Loop checking, 11.06.2015, T. Koi
80 
81  return( ans );
82 }
const int N
Definition: mixmax.h:43
const char * p
Definition: xmltok.h:285

Here is the caller graph for this function:

double nf_polevl ( double  x,
double  coef[],
int  N 
)

Definition at line 46 of file nf_polevl.cc.

46  {
47 
48  double ans;
49  int i;
50  double *p;
51 
52  p = coef;
53  ans = *p++;
54  i = N;
55 
56  do {
57  ans = ans * x + *p++; }
58  while( --i ); // Loop checking, 11.06.2015, T. Koi
59 
60  return( ans );
61 }
const int N
Definition: mixmax.h:43
const char * p
Definition: xmltok.h:285

Here is the caller graph for this function: