Main Page | Namespace List | Class Hierarchy | Class List | File List | Class Members | File Members | Related Pages

math_util.cpp

Go to the documentation of this file.
00001 /*
00002  * ******** fete: From ENDF To ENDL *********
00003  * 
00004  * Copyright (c) 2006, The Regents of the University of California. 
00005  * All rights reserved.
00006  * 
00007  * Produced at the Lawrence Livermore National Laboratory. 
00008  * Written by David A. Brown, Gerry Hedstrom, Tony Hill
00009  * 
00010  * This file is part of fete v1.0  (UCRL-CODE-218718)
00011  * 
00012  * Please read the COPYING file for "Our Notice and GNU General 
00013  * Public License" in the root of this software distribution.  
00014  * 
00015  * This program is free software; you can redistribute it and/or modify 
00016  * it under the terms of the GNU General Public License (as published by 
00017  * the Free Software Foundation) version 2, dated June 1991. 
00018  * 
00019  * This program is distributed in the hope that it will be useful, 
00020  * but WITHOUT ANY WARRANTY; without even the IMPLIED WARRANTY OF 
00021  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the terms 
00022  * and conditions of the GNU General Public License for more details. 
00023  * 
00024  * You should have received a copy of the GNU General Public License along 
00025  * with this program; if not, write to the Free Software Foundation, Inc., 
00026  * 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 
00027  * 
00028  * $Revision: 1872 $
00029  * $Date: 2006-05-17 08:34:30 -0700 (Wed, 17 May 2006) $
00030  * $Author: dbrown $
00031  * $Id: math_util.cpp 1872 2006-05-17 15:34:30Z dbrown $
00032  * 
00033  * ******** fete: From ENDF To ENDL *********
00034  */
00035 
00036 // implementation of mymatrix class
00037 
00038 #include <stdio.h>
00039 #include <stdlib.h>
00040 #include <math.h>
00041 #include <float.h>
00042 
00043 #include "math_util.hpp"
00044 #include "messaging.hpp"
00045 #include "ENDF_file.hpp"
00046 
00047 using namespace std;
00048 extern ENDLClass ENDL;
00049 
00050 // *************** Param class *******************
00051 
00052 // ------------- Param::set_zero -------------
00053 void Param::set_zero( )
00054 // Set parameters to zero
00055 {
00056   for(Param::iterator itr=begin(); itr!=end(); ++itr)
00057   {
00058     *itr = 0.0;
00059   }
00060 }
00061 
00062 // ------------- Param::print() -------------
00063 void Param::print()
00064 {
00065   cout.setf(ios::scientific,ios::floatfield);
00066   for(Param::iterator itr=begin(); itr!=end(); ++itr)
00067   {
00068     cout << ENDL.data( *itr ) << "  ";
00069   }
00070   cout << endl;
00071 }
00072 
00073 // ------------- Param::operator*= -------------
00074 Param& Param::operator*=( double alpha )
00075 {
00076   for(Param::iterator itr=begin(); itr!=end(); ++itr)
00077   {
00078     *itr *= alpha;
00079   }
00080   return *this;
00081 }
00082 
00083 // ------------- Param::operator+= -------------
00084 Param& Param::operator+=( Param& to_add )
00085 {
00086   if( size() != to_add.size() )
00087   {
00088     SevereError( "Param::operator+=", "incompatible sizes" );
00089   }
00090   for( unsigned int count = 0; count < size(); ++count )
00091   {
00092     at( count ) += to_add[ count ];
00093   }
00094   return *this;
00095 }
00096 
00097 // *************** basic interpolation routines *****************
00098 // ------------------- Lin_Lin ------------------------
00099 double Lin_Lin( double x_mid, double x_0, double y_0,
00100   double x_1, double y_1 )
00101 {
00102   double alpha = ( x_mid - x_0 ) / ( x_1 - x_0 );
00103   if((alpha <= 0.0) || (alpha >= 1.0))
00104   {
00105     SevereError("Lin_Lin","Attempt to extrapolate");
00106   }
00107   return (1.0 - alpha)*y_0 + alpha*y_1;
00108 }
00109 // ------------------- Log_Lin ------------------------
00110 double Log_Lin( double x_mid, double x_0, double y_0,
00111   double x_1, double y_1 )
00112 {
00113   double x_ratio = x_1 / x_0;
00114 
00115   // check for division by 0 and bad logarithm
00116   if( (x_ratio <= 0.0) || (x_ratio == 1.0) )
00117   {
00118     SevereError("Log_Lin",pastenum("bad x_ratio:",x_ratio));
00119   }
00120 
00121   double b = ( y_1 - y_0 ) / log( x_ratio );
00122 
00123   // check for extrapolation
00124   if( ( (x_mid > x_0) && (x_mid > x_1) ) ||
00125       ( (x_mid < x_0) && (x_mid < x_1) ) )
00126   {
00127     SevereError("Log_Lin","Attempt to extrapolate");
00128   }
00129   return y_0 + b*log( x_mid / x_0 );
00130 }
00131 // ------------------- Lin_Log ------------------------
00132 double Lin_Log( double x_mid, double x_0, double y_0,
00133   double x_1, double y_1 )
00134 // interpolation for log(y) = log(y_0) + b*(x - x_0)
00135 {
00136   // check for proper inputs
00137   if( (y_0 <= 0.0) || (y_1 <= 0.0) ||
00138       (x_1 == x_0) )
00139   {
00140     if ( ( x_mid == x_0 ) && ( y_0 == y_1 ) ) return y_0;
00141     SevereError("Lin_Log", 
00142         "Bad inputs, (x_0,y_0), (x_1,y_1), x_mid: "+
00143         pastenum("(",x_0)+
00144         pastenum(",",y_0)+
00145         pastenum("), (",x_1)+
00146         pastenum(",",y_1)+
00147         pastenum("), ",x_mid)
00148     );
00149   }
00150   double b = log( y_1 / y_0 ) / ( x_1 - x_0 );
00151 
00152   // check for extrapolation
00153   if( ( (x_mid > x_0) && (x_mid > x_1) ) ||
00154       ( (x_mid < x_0) && (x_mid < x_1) ) )
00155   {
00156     SevereError("Lin_Log","Attempt to extrapolate");
00157   }
00158 
00159   return ( y_0 * exp( b * ( x_mid - x_0 ) ) );
00160 }
00161 // ------------------- Log_Log ------------------------
00162 double Log_Log( double x_mid, double x_0, double y_0,
00163   double x_1, double y_1 )
00164 // interpolation for log(y) = log(y_0) + b*log(x/x_0)
00165 {
00166   if( (y_0 <= 0.0) || (y_1 <= 0.0) ||
00167       (x_0 <= 0.0) || (x_1 <= 0.0) ||
00168       (x_1 == x_0) )
00169   {
00170     if ( ( x_mid == x_0 ) && ( y_0 == y_1 ) && ( x_0 == x_1 ) ) return y_0;
00171     SevereError("Log_Log","Bad inputs");
00172   }
00173 
00174   double b = log( y_1 / y_0 ) /
00175     log( x_1 / x_0 );
00176 
00177   // check for extrapolation
00178   if( ( (x_mid > x_0) && (x_mid > x_1) ) ||
00179       ( (x_mid < x_0) && (x_mid < x_1) ) )
00180   {
00181     SevereError("Log_Log","Attempt to extrapolate");
00182   }
00183 
00184   return y_0 * exp( b * log(x_mid/x_0) );
00185 }
00186 
00187 // ---------------- maxE_loglog ----------------------------
00188 double maxE_loglog( double xi, double alpha )
00189 // get the maximum error for linear interpolation of x^alpha
00190 // on the interval 1 < x < xi.
00191 // This routine is used in linear interpolation of loglog data.
00192 {
00193   if ( xi <= 1.0) //Make sure xi is not below interval
00194   {
00195     SevereError("maxE_loglog","xi must be larger than 1!");
00196   }
00197 
00198   //Calculate x_tilde - place where error is largest in interval
00199   double x_tilde;
00200   double u_tilde;  // y value at x_tilde
00201   double log_xi = log(xi);
00202 
00203   if ( abs(alpha*log_xi) <= 0.00001 )  //Take care of pole at 0
00204   {
00205     x_tilde = (1.0 - xi + alpha*log_xi)/((alpha - 1.0)*log_xi);
00206   }
00207   else if ( abs((alpha-1.0)*log_xi) <= 0.00001 )  //Take care of pole at +1
00208   {
00209     x_tilde = (alpha*xi*log_xi)/(xi - 1.0 + (alpha-1.0)*xi*log_xi);
00210   }
00211   else
00212   {
00213     x_tilde =
00214       alpha*(pow(xi,alpha) - xi)/((alpha - 1.0)*(pow(xi,alpha) - 1.0));
00215   }
00216   u_tilde = pow(x_tilde, alpha);
00217   //  cerr << "x_tilde = " << x_tilde << "  u_tilde = " << u_tilde << endl;
00218 
00219   //Comput v_tilde (and beta) --- y value of linlin extrap at x_tilde
00220   double beta = (pow(xi,alpha)-1.0)/(xi-1.0);
00221   //  cerr<<"beta = "<<beta<<endl;
00222   double v_tilde = 1.0 + beta*(x_tilde - 1.0);
00223 
00224   //Relative error at x_tilde
00225   double rel_err = abs((v_tilde - u_tilde)/u_tilde);
00226   //  cerr << "v_tilde = " << v_tilde << "  rel_err = " << rel_err << endl;
00227 
00228   return rel_err;
00229 }
00230 
00231 // ******* translation of Slatec incomplete gamma function *****
00232 // local variables:
00233 static int KONTRL = 2;
00234 static char StrLevel[4][32] = { { "Informative message" },
00235   { "Informative message" },
00236   { "recoverable error" },
00237   { "fatal error" } };
00238 
00239 // local routines
00240 double c_d1mach( int I );
00241 double c_d9gmit( double a, double x, double algap1, double sgngam, double alx );
00242 double c_d9lgic( double a, double x, double alx );
00243 double c_d9lgit( double a, double x, double algap1 );
00244 double c_d9lgmc( double x );
00245 double c_dcsevl( double x, double *cs, int n );
00246 void c_dgamlm( double *xmin, double *xmax );
00247 double c_dgamr( double x );
00248 void c_dlgams( double x, double *dlgam, double *sgngam );
00249 double c_dlngam( double x );
00250 int c_initds( double *os, int nos, double eta );
00251 void c_xermsg( char *Lib, char *Routine, char *Msg, int nError, int level );
00252 void c_xsetf( int k ); 
00253 int c_xgetf( void );
00254 static void c_xerrPrint( char *Lib, char *Routine, char *Msg, int nError,
00255   int level );
00256 
00257 /*
00258 ********************** c_d1mach *****************************
00259 */
00260 double c_d1mach( int I )
00261 // machine-dependent parameters
00262 {
00263 
00264         static double D1MACH[5] = { DBL_MIN, DBL_MAX, DBL_EPSILON / 2.,
00265                            DBL_EPSILON, 3.01029995663981198e-01 };
00266 
00267         if( ( I < 1 ) || ( I > 5 ) ) c_xermsg( "cSLATEC", "d1mach",
00268            "I out of bounds", 1, 2 );
00269 
00270         return( D1MACH[I-1] );
00271 }
00272 
00273 /*
00274 *********************** c_d9gmit ********************************
00275 */
00276 double c_d9gmit( double a, double x, double algap1, double sgngam,
00277   double alx ) 
00278 // Compute Tricomi's incomplete Gamma function for small arguments.
00279 {
00280 
00281         static int FirstCall = 1;
00282         static double eps, bot;
00283         int k, m, ma;
00284         double ae, aeps, algs, alg2, fk, s, sgng2, t, te, d9gmit;
00285 
00286         if( FirstCall ) {
00287                 FirstCall = 0;
00288                 eps = 0.5 * c_d1mach( 3 );
00289                 bot = log( c_d1mach( 1 ) );
00290         }
00291 
00292         if( x <= 0. ) c_xermsg( "cSLATEC", "c_d9gmit", "X should be > 0", 1, 2 );
00293 
00294         ma = static_cast<int>( a + 0.5 );
00295         if( a < 0. ) ma--;
00296         aeps = a - ma;
00297 
00298         ae = a;
00299         if( a < -0.5 ) ae = aeps;
00300 
00301         t = 1.;
00302         te = ae;
00303         s = t;
00304         for( k = 1; k <= 200; k++ ) {
00305                 fk = k;
00306                 te = -x * te / fk;
00307                 t = te / ( ae + fk );
00308                 s = s + t;
00309                 if( fabs( t ) < eps * fabs( s ) ) break;
00310         }
00311         if( k > 200 ) c_xermsg( "cSLATEC", "c_d9gmit",
00312          "No convergence in 200 terms of taylor-s series", 2, 2 );
00313 
00314         if( a >= -0.5 ) {
00315                 algs = -algap1 + log(s);
00316                 d9gmit = exp( algs ); }
00317         else {
00318                 algs = -c_dlngam( 1. + aeps ) + log( s );
00319                 s = 1.;
00320                 m = -ma - 1;
00321                 if( m != 0 ) {
00322                         t = 1.;
00323                         for( k = 0; k < m; k++ ) {
00324                                 t = x * t / ( aeps - ( m - k ) );
00325                                 s = s + t;
00326                                 if( fabs( t ) < eps * fabs( s ) ) break;
00327                         }
00328                 }
00329                 d9gmit = 0.;
00330                 algs = -ma * log( x ) + algs;
00331                 if( ( s == 0. ) && ( aeps == 0. ) ) return( exp( algs ) );
00332                 sgng2 = 1.;
00333                 if( s < 0. ) sgng2 = -1.;
00334                 alg2 = -x - algap1 + log( fabs ( s ) );
00335                 if( alg2 > bot ) d9gmit = sgng2 * exp( alg2 );
00336                 if( algs > bot ) d9gmit = d9gmit + exp( algs );
00337         }
00338         return( d9gmit );
00339 }
00340 
00341 /*
00342 ****************************** c_d9lgic ****************************
00343 */
00344 double c_d9lgic( double a, double x, double alx ) 
00345 // Compute the log complementary incomplete Gamma function
00346 //            for large X and for A .LE. X.
00347 {
00348 
00349         static double eps = DBL_EPSILON / 2.;
00350         int k;
00351         double fk, p, r, s, t, xma, xpa;
00352 
00353         xpa = x + 1. - a;
00354         xma = x - 1. - a;
00355 
00356         r = 0.;
00357         p = 1.;
00358         s = p;
00359         for( k = 1; k <= 300; k++ ) {
00360                 fk = k;
00361                 t = fk * ( a - fk ) * ( 1. + r );
00362                 r = -t / ( ( xma + 2. * fk ) * ( xpa + 2. * fk ) + t );
00363                 p = r * p;
00364                 s = s + p;
00365                 if( fabs( p ) < eps * s ) break;
00366         }
00367         if( k > 300 ) c_xermsg ( "cSLATEC", "c_d9lgic", 
00368             "no convergence in 300 terms of continued fraction", 1, 2 );
00369         return( a * alx - x + log( s / xpa ) );
00370 }
00371 
00372 /*
00373 ****************************** c_d9lgit ******************************
00374 */
00375 double c_d9lgit( double a, double x, double algap1 ) 
00376 // Compute the logarithm of Tricomi's incomplete Gamma
00377 //            function with Perron's continued fraction for large X and
00378 //            A .GE. X.
00379 {
00380 
00381         static int FirstCall = 1;
00382         static double eps, sqeps;
00383         int k;
00384         double ax, a1x, fk, hstar, p, r, s, t;
00385 
00386         if( FirstCall ) {
00387                 FirstCall = 0;
00388                 eps = 0.5 * c_d1mach( 3 );
00389                 sqeps = sqrt( c_d1mach( 4 ) );
00390         }
00391 
00392         if( ( x <= 0. ) || ( a < x ) ) c_xermsg( "cSLATEC", "c_d9lgit",
00393            "x should be > 0.0 and <= a", 2, 2 );
00394 
00395         ax = a + x;
00396         a1x = ax + 1.0;
00397         r = 0.;
00398         p = 1.;
00399         s = p;
00400         for( k = 1; k <= 200; k++ ) {
00401                 fk = k;
00402                 t = ( a + fk ) * x * ( 1. + r );
00403                 r = t / ( ( ax + fk ) * ( a1x + fk ) - t );
00404                 p = r * p;
00405                 s = s + p;
00406                 if( fabs( p ) < eps * s) break;
00407         }
00408         if( k > 200 ) c_xermsg ( "cSLATEC", "c_d9lgit",
00409              "no convergence in 200 terms of continued fraction", 3, 2 );
00410         hstar = 1.0 - x * s / a1x;
00411         if( hstar < sqeps) c_xermsg ( "cSLATEC", "c_d9lgit",
00412                "result less than half precision", 1, 1 );
00413         return( -x - algap1 - log( hstar ) );
00414 }
00415 
00416 /*
00417 *************************** c_d9lgmc *********************************
00418 */
00419 double c_d9lgmc( double x ) 
00420 // Compute the log Gamma correction factor so that
00421 //            LOG(DGAMMA(X)) = LOG(SQRT(2*PI)) + (X-5.)*LOG(X) - X
00422 //            + D9LGMC(X).
00423 {
00424 
00425 
00426         static int FirstCall = 1;
00427         static int nalgm;
00428         static double xbig, xmax, algmcs[15] = {
00429                 +0.1666389480451863247205729650822e+00,
00430                 -0.1384948176067563840732986059135e-04,
00431                 +0.9810825646924729426157171547487e-08,
00432                 -0.1809129475572494194263306266719e-10,
00433                 +0.6221098041892605227126015543416e-13,
00434                 -0.3399615005417721944303330599666e-15,
00435                 +0.2683181998482698748957538846666e-17,
00436                 -0.2868042435334643284144622399999e-19,
00437                 +0.3962837061046434803679306666666e-21,
00438                 -0.6831888753985766870111999999999e-23,
00439                 +0.1429227355942498147573333333333e-24,
00440                 -0.3547598158101070547199999999999e-26,
00441                 +0.1025680058010470912000000000000e-27,
00442                 -0.3401102254316748799999999999999e-29,
00443                 +0.1276642195630062933333333333333e-30 };
00444         double d9lgmc;
00445 
00446         if( FirstCall ) {
00447                 FirstCall = 0;
00448                 nalgm = c_initds( algmcs, 15, c_d1mach( 3 ) );
00449                 xbig = 1.0 / sqrt( c_d1mach( 3 ) );
00450                 xmax = c_d1mach( 2 ) / 12.;
00451                 if( xmax > 1. / ( 12. * c_d1mach( 1 ) ) )
00452                    xmax = 1. / ( 12. * c_d1mach( 1 ) );
00453                 xmax = exp( xmax );
00454         }
00455 
00456         if( x < 10. ) c_xermsg( "cSLATEC", "c_d9lgmc", "x must be >= 10", 1, 2 );
00457         if( x < xmax ) {
00458                 d9lgmc = 1. / ( 12. * x );
00459                 if( x < xbig ) d9lgmc =
00460                   c_dcsevl( 2.0 * ( 10. / x ) * ( 10. / x ) - 1., algmcs, nalgm ) / x;
00461         }
00462         else {
00463                 d9lgmc = 0.;
00464                 c_xermsg( "cSLATEC", "c_d9lgmc",
00465                     "x so big d9lgmc underflows", 2, 1 );
00466         }
00467         return( d9lgmc );
00468 }
00469 
00470 /*
00471 ************************* c_dcsevl ***********************************
00472 */
00473 double c_dcsevl( double x, double *cs, int n ) 
00474 // Evaluate the N-term Chebyshev series CS at X.
00475 {
00476 
00477         static int FirstCall = 1;
00478         static double onepl;
00479         int i, ni;
00480         double b0, b1, b2, twox;
00481 
00482         if( FirstCall ) {
00483                 FirstCall = 0;
00484                 onepl = 1.0 + c_d1mach( 4 );
00485         }
00486 
00487         if( n < 1 ) c_xermsg( "cSLATEC", "c_dcsevl",
00488              "number of terms <= 0", 2, 2 );
00489         if( n > 1000 ) c_xermsg( "cSLATEC", "c_dcsevl",
00490              "number of terms > 1000", 3, 2 );
00491         if( fabs( x ) > onepl ) c_xermsg( "cSLATEC", "c_dcsevl",
00492              "x outside the interval (-1, +1)", 1, 1 );
00493 
00494         b1 = 0.;
00495         b0 = 0.;
00496         twox = 2.0 * x;
00497         for( i = 0; i <= n; i++ ) {
00498                 b2 = b1;
00499                 b1 = b0;
00500                 ni = n - i;
00501                 b0 = twox * b1 - b2 + cs[ ni ];
00502         }
00503         return( 0.5 * ( b0 - b2 ) );
00504 }
00505 
00506 /*
00507 ************************** c_dgamit **********************************
00508 */
00509 double c_dgamit( double a, double x ) 
00510 // Calculate Tricomi's form of the incomplete Gamma function,
00511 //   x^{-a} / Gamma(a) * int_0^x t^(a-1) * exp{-t} dt.
00512 {
00513 
00514         static int FirstCall = 1;
00515         static double alneps, sqeps;//, bot;
00516         double aeps, ainta, algap1, alng, alx, h, sga, sgngam, t, dgamit;
00517 
00518         if( FirstCall ) {
00519                 FirstCall = 0;
00520                 alneps = -log( c_d1mach( 3 ) );
00521                 sqeps = sqrt( c_d1mach( 4 ) );
00522                 //bot = log( c_d1mach( 1 ) );
00523         }
00524 
00525         if( x < 0. ) c_xermsg( "cSLATEC", "c_dgamit", "x is negative", 2, 2 );
00526         if( x != 0. ) alx = log( x );
00527         sga = 1.;
00528         if( a < 0. ) sga = -1.;
00529         ainta = (double) ( (int) ( a + 0.5 * sga ) );
00530         aeps = a - ainta;
00531 
00532         if( x <= 0. ) {
00533                 dgamit = 0.;
00534                 if( ( ainta > 0. ) || ( aeps != 0. ) )
00535                     dgamit = c_dgamr( a + 1.0 ); }
00536         else if( x <= 1. ) {
00537                         if( ( a >= -0.5 ) || ( aeps != 0. ) )
00538                             c_dlgams( a + 1., &algap1, &sgngam );
00539                         dgamit = c_d9gmit( a, x, algap1, sgngam, alx ); }
00540         else if( a >= x ) {
00541                 t = c_d9lgit( a, x, c_dlngam( a + 1. ) );
00542                 dgamit = exp( t ); }
00543         else {
00544                 alng = c_d9lgic( a, x, alx );
00545                 h = 1.;
00546                 if( ( aeps != 0. ) || ( ainta > 0. ) ) {
00547                         c_dlgams( a + 1., &algap1, &sgngam );
00548                         t = log( fabs( a ) ) + alng - algap1;
00549                         if( t > alneps ) goto l60;
00550                         if( t > -alneps) h = 1. - sga * sgngam * exp( t );
00551                         if( fabs( h ) <= sqeps )
00552                            c_xermsg( "cSLATEC", "c_dgamit",
00553                                     "result lt half precision", 1, 1 );
00554                 }
00555                 t = -a * alx + log( fabs( h ) );
00556                 dgamit = exp( t );
00557                 if( h < 0. ) dgamit = -dgamit;
00558                 return( dgamit );
00559 
00560 l60:    t = t - a * alx;
00561                 dgamit = -sga * sgngam * exp( t );
00562         }
00563         return( dgamit );
00564 }
00565 
00566 /*
00567 *********************** c_dgamlm *************************************
00568 */
00569 void c_dgamlm( double *xmin, double *xmax ) 
00570 // Compute the minimum and maximum bounds for the argument in
00571 //            the Gamma function.
00572 {
00573 
00574         int i;
00575         double alnbig, alnsml, xln, xold;
00576 
00577         alnsml = log( c_d1mach( 1 ) );
00578         *xmin = -alnsml;
00579         for( i = 0; i < 10; i++ ) {
00580                 xold = *xmin;
00581                 xln = log( *xmin );
00582                 *xmin = *xmin - *xmin * ( ( *xmin + 0.5 ) * xln - *xmin - 0.2258
00583  + alnsml ) / ( *xmin * xln + 0.5 );
00584                 if( fabs( *xmin - xold ) < 0.005 ) break;
00585         }
00586         if( i == 10 )
00587            c_xermsg( "cSLATEC", "c_dgamlm", "unable to find xmin", 1, 2 );
00588         *xmin = -*xmin + 0.01;
00589         alnbig = log( c_d1mach( 2 ) );
00590         *xmax = alnbig;
00591         for( i = 0; i < 10; i++ ) {
00592                 xold = *xmax;
00593                 xln = log( *xmax );
00594                 *xmax = *xmax - *xmax * ( ( *xmax - 0.5 ) * xln - *xmax + 0.9189
00595  - alnbig ) / ( *xmax * xln - 0.5 );
00596                 if( fabs( *xmax - xold ) < 0.005 ) break;
00597         }
00598         if( i == 10 )
00599            c_xermsg( "cSLATEC", "c_dgamlm", "unable to find xmax", 2, 2 );
00600         *xmax = *xmax - 0.01;
00601         if( ( -*xmax + 1. ) > *xmin ) *xmin = -*xmax + 1.;
00602 }
00603 /*
00604 **************************** c_dgamma ********************************
00605 */
00606 double c_dgamma( double x ) 
00607 // Compute the complete Gamma function.
00608 {
00609 
00610         static int FirstCall = 1, ngam;
00611         static double pi = 3.14159265358979323846264338327950,
00612                 sq2pil = 0.91893853320467274178032973640562,
00613                  xmin, xmax, dxrel;
00614         static double gamcs[42] = {
00615                 +0.8571195590989331421920062399942e-02,
00616                 +0.4415381324841006757191315771652e-02,
00617                 +0.5685043681599363378632664588789e-01,
00618                 -0.4219835396418560501012500186624e-02,
00619                 +0.1326808181212460220584006796352e-02,
00620                 -0.1893024529798880432523947023886e-03,
00621                 +0.3606925327441245256578082217225e-04,
00622                 -0.6056761904460864218485548290365e-05,
00623                 +0.1055829546302283344731823509093e-05,
00624                 -0.1811967365542384048291855891166e-06,
00625                 +0.3117724964715322277790254593169e-07,
00626                 -0.5354219639019687140874081024347e-08,
00627                 +0.9193275519859588946887786825940e-09,
00628                 -0.1577941280288339761767423273953e-09,
00629                 +0.2707980622934954543266540433089e-10,
00630                 -0.4646818653825730144081661058933e-11,
00631                 +0.7973350192007419656460767175359e-12,
00632                 -0.1368078209830916025799499172309e-12,
00633                 +0.2347319486563800657233471771688e-13,
00634                 -0.4027432614949066932766570534699e-14,
00635                 +0.6910051747372100912138336975257e-15,
00636                 -0.1185584500221992907052387126192e-15,
00637                 +0.2034148542496373955201026051932e-16,
00638                 -0.3490054341717405849274012949108e-17,
00639                 +0.5987993856485305567135051066026e-18,
00640                 -0.1027378057872228074490069778431e-18,
00641                 +0.1762702816060529824942759660748e-19,
00642                 -0.3024320653735306260958772112042e-20,
00643                 +0.5188914660218397839717833550506e-21,
00644                 -0.8902770842456576692449251601066e-22,
00645                 +0.1527474068493342602274596891306e-22,
00646                 -0.2620731256187362900257328332799e-23,
00647                 +0.4496464047830538670331046570666e-24,
00648                 -0.7714712731336877911703901525333e-25,
00649                 +0.1323635453126044036486572714666e-25,
00650                 -0.2270999412942928816702313813333e-26,
00651                 +0.3896418998003991449320816639999e-27,
00652                 -0.6685198115125953327792127999999e-28,
00653                 +0.1146998663140024384347613866666e-28,
00654                 -0.1967938586345134677295103999999e-29,
00655                 +0.3376448816585338090334890666666e-30,
00656                 -0.5793070335782135784625493333333e-31 };
00657         int i, n;
00658         double sinpiy, y, dgamma;
00659 
00660         if( FirstCall ) {
00661                 FirstCall = 0;
00662                 ngam = c_initds( gamcs, 42, 0.1 * c_d1mach( 3 ) );
00663                 c_dgamlm( &xmin, &xmax );
00664                 dxrel = sqrt( c_d1mach( 4 ) );
00665         }
00666 
00667         y = fabs( x );
00668         if( y <= 10 ) {
00669                 n = static_cast<int>( x );
00670                 if( x < 0. ) n--;
00671                 y = x - n;
00672                 n--;
00673                 dgamma = 0.9375 + c_dcsevl( 2. * y - 1., gamcs, ngam );
00674                 if( n != 0 ) {
00675                         if( n < 0 ) {
00676                                 n = -n;
00677                                 if( x == 0. )
00678                                    c_xermsg( "cSLATEC", "c_dgamma",
00679                                        "x is 0", 4, 2 );
00680                                 if( ( x <  0.0 ) && ( x + n - 2 == 0. ) )
00681                                     c_xermsg( "cSLATEC", "c_dgamma", 
00682                                       "x is a negative integer", 4, 2 );
00683                                 if( ( x < -0.5 ) && ( fabs( ( x - ( (int) ( x - 
00684 0.5 ) ) ) / x ) < dxrel ) )
00685                                      c_xermsg( "cSLATEC", "c_dgamma",
00686                 "answer < half precision because x too near negative integer", 1, 1 );
00687                                 for( i = 0; i < n; i++ ) dgamma /= ( x + i ); }
00688                         else {
00689                                 for( i = 1; i <= n; i++ ) dgamma *= ( y + i );
00690                         }
00691                 } }
00692         else {
00693                 if( x > xmax )
00694                     c_xermsg( "cSLATEC", "c_dgamma",
00695                          "x so big gamma overflows", 3, 2 );
00696                 dgamma = 0.;
00697                 if(x  < xmin )
00698                     c_xermsg( "cSLATEC", "c_dgamma", "u so small gamma underflows",
00699                        2, 1 );
00700                 if(x >= xmin ) {
00701                         dgamma = exp( ( y - 0.5 ) * log( y ) - y + sq2pil + c_d9lgmc( y ) );
00702                         if( x <= 0. ) {
00703                                 if( fabs( ( x - (int) (x - 0.5 ) ) / x ) < dxrel
00704  )
00705                                  c_xermsg( "cSLATEC", "c_dgammae",
00706                     "answer lt half precision, x too near negative integer", 1, 1 );
00707                                 sinpiy = sin( pi * y );
00708                                 if( sinpiy == 0. )
00709                       c_xermsg( "c_SLATEC", "c_dgamma", "x is a negative integer",
00710                             4, 2 );
00711                                 dgamma = -pi / ( y * sinpiy * dgamma );
00712                         }
00713                 }
00714         }
00715         return( dgamma );
00716 }
00717 
00718 /*
00719 ************************* c_dgamr ***********************************
00720 */
00721 double c_dgamr( double x ) 
00722 // Compute the reciprocal of the Gamma function.
00723 {
00724 
00725         int irold;
00726         double dgamr, alngx, sgngx;
00727 
00728         dgamr = 0.;
00729         if( ( x > 0. ) || ( (double) ( (int) x ) != x ) ) {
00730                 irold = c_xgetf( );
00731                 c_xsetf( 1 );
00732                 if( fabs( x ) <= 10. ) {
00733                         dgamr = 1. / c_dgamma( x );
00734                         c_xsetf( irold ); }
00735                 else {
00736                         c_dlgams( x, &alngx, &sgngx );
00737                         c_xsetf( irold );
00738                         dgamr = sgngx * exp( -alngx );
00739                 }
00740         }
00741         return( dgamr );
00742 }
00743 
00744 /*
00745 *************************** c_dlgams *********************************
00746 */
00747 void c_dlgams( double x, double *dlgam, double *sgngam ) 
00748 // Compute the logarithm of the absolute value of the Gamma
00749 //            function.
00750 {
00751 
00752         *dlgam = c_dlngam( x );
00753         *sgngam = 1.;
00754         if( x <= 0. ) {
00755                 if( ( ( -(int) x ) % 2 ) == 0 ) *sgngam = -1.;
00756         }
00757 }
00758 
00759 /*
00760 ************************ c_dlngam ************************************
00761 */
00762 double c_dlngam( double x ) 
00763 // Compute the logarithm of the absolute value of the Gamma
00764 //            function.
00765 {
00766 
00767         static int FirstCall = 1;
00768         static double sq2pil = 0.91893853320467274178032973640562,
00769                   sqpi2l = 0.225791352644727432363097614947441,
00770                   pi     = 3.14159265358979323846264338327950, xmax, dxrel;
00771         double sinpiy, y, temp;
00772 
00773         if( FirstCall ) {
00774                 FirstCall = 0;
00775                 temp = 1. / log( c_d1mach( 2 ) );
00776                 xmax = temp * c_d1mach( 2 );
00777                 dxrel = sqrt( c_d1mach( 4 ) );
00778         }
00779 
00780         y = fabs( x );
00781         if( y <= 10. ) return( log( fabs( c_dgamma( x ) ) ) );
00782 
00783         if( y > xmax )
00784              c_xermsg( "cSLATEC", "c_dlngam", "abs(x) so big c_dlngam overflows",
00785                    2, 2 );
00786         if( x > 0 ) return( sq2pil + ( x - 0.5 ) * log( x ) - x + c_d9lgmc( y ) 
00787 );
00788 
00789         sinpiy = fabs( sin( pi * y ) );
00790         if( sinpiy == 0. )
00791               c_xermsg( "cSLATEC", "c_dlngam", "x is a negative integer", 3, 2 );
00792         if( fabs( ( x - (int) ( x -0.5 ) ) / x ) < dxrel )
00793               c_xermsg( "c_SLATEC", "c_dlngam",
00794                 "answer < half precision because x too near negative integer",
00795                      1, 1 );
00796         return( sqpi2l + ( x - 0.5 ) * log( y ) - x - log( sinpiy ) - c_d9lgmc( y ) );
00797 }
00798 
00799 /*
00800 **************************** c_initds ********************************
00801 */
00802 int c_initds( double *os, int nos, double eta ) 
00803 // Determine the number of terms needed in an orthogonal
00804 //            polynomial series so that it meets a specified accuracy.
00805 {
00806 
00807         int i, ii;
00808         double err;
00809 
00810         if( nos < 1 )
00811             c_xermsg( "cSLATEC", "c_initds", "number of coefficients is less than 1",
00812                 2, 1 );
00813         err = 0.;
00814         for( ii = 1; ii <= nos; ii++ ) {
00815                 i = nos - ii;
00816                 err += fabs( os[i] );
00817                 if( err > eta ) break;
00818         }
00819         if( i >= nos )
00820             c_xermsg( "cSLATEC", "c_initds",
00821                 "Chebyshev series too short for specified accuracy", 1, 1 );
00822         return( i );
00823 } 
00824 
00825 /*
00826 ************************ c_xermsg **********************************
00827 */
00828 void c_xermsg( char *Lib, char *Routine, char *Msg, int nError, int level ) 
00829 // Process an error message depending on the severity level
00830 {
00831 
00832         char Str[1024];
00833 
00834         if( ( level < -1 ) || ( level > 2 ) ) {
00835                 sprintf( Str, "invalid level = %d", level );
00836                 c_xerrPrint( Lib, "c_xermsg", Str, 1, 2 );
00837                 exit( EXIT_FAILURE );
00838         }
00839 
00840         c_xerrPrint( Lib, Routine, Msg, nError, level );
00841         if( level == 2 ) {
00842                 exit( EXIT_FAILURE );
00843         }
00844 }
00845 /*
00846 ************************** c_xerrPrint ********************************
00847 */
00848 void c_xerrPrint( char *Lib, char *Routine, char *Msg, int nError, int level ) 
00849 // Print an error message
00850 {
00851 
00852         fprintf( stderr, "\n" );
00853         fprintf( stderr, "*** Message from routine %s in library %s\n", Routine,
00854  Lib );
00855         fprintf( stderr, "*** %s\n", StrLevel[level+1] );
00856         fprintf( stderr, "*   %s\n", Msg );
00857         fprintf( stderr, "*   Error number = %d\n", nError );
00858         fprintf( stderr, "*\n" );
00859         fprintf( stderr, "*** End of message\n" );
00860         fprintf( stderr, "\n" );
00861 }
00862 
00863 /*
00864 *********************** c_xsetf ***********************************
00865 */
00866 void c_xsetf( int k ) 
00867 // set the error severity level to k
00868 {
00869 
00870         char Str[1024];
00871 
00872         if( abs( k ) > 2 ) {
00873                 sprintf( Str, "invalid k = %d value", k );
00874                 c_xerrPrint( "cSlatec", "c_xsetf", Str, 1, 2 );
00875         }
00876         KONTRL = abs( k );
00877         return;
00878 }
00879 
00880 /*
00881 ****************************** c_xgetf **************************
00882 */
00883 int c_xgetf( void ) { return( KONTRL ); }
00884 
00885 // ************* expint function *******************************
00886 /*********************************************************************
00887    Returns the exponential integral function
00888    E_n(x) = int_1^infinity e^{-x*t}/t^n dt,     for x > 0.
00889    C.A. Bertulani        May/15/2000
00890 *********************************************************************/
00891 #define EULER 0.5772156649      /* Euler's constant gamma */
00892 #define MAXIT 100               /* Maximum allowed number of iterations. */
00893 #define FPMIN 1.0e-30  // close to the smallest representable floting-point number.
00894 #define EPS 6.0e-8    // Desired relative error, not smaller than the machine precision.
00895 
00896 double expint(int n, double x)
00897 {
00898   int i,ii,nm1;
00899   double a,b,c,d,del,fact,h,psi,ans;
00900 
00901   nm1=n-1;
00902   if (n < 0 || x < 0.0 || (x==0.0 && (n==0 || n==1)))
00903   {
00904     SevereError("expint", "Bad arguments");
00905   }
00906   else
00907   {
00908     if (n == 0)
00909     {
00910       ans=exp(-x)/x;   /* Special case */
00911     }
00912     else if (x == 0.0)
00913     {
00914       ans=1.0/nm1;  /* Another special case */
00915     }
00916     else if (x > 1.0)
00917     {   /* Lentz's algorithm */
00918       b=x+n;
00919       c=1.0/FPMIN;
00920       d=1.0/b;
00921       h=d;
00922       for (i=1;i<=MAXIT;i++)
00923       {
00924         a = -i*(nm1+i);
00925         b += 2.0;
00926         d=1.0/(a*d+b);  /* Denominators cannot be zero */
00927         c=b+a/c;
00928         del=c*d;
00929         h *= del;
00930         if (fabs(del-1.0) < EPS)
00931         {
00932           ans=h*exp(-x);
00933           return ans;
00934         }
00935       }
00936       SevereError("expint", "Continued fraction failed");
00937     }
00938     else
00939     {
00940       ans = (nm1!=0) ? 1.0/nm1 : -log(x)-EULER; /* Set first term */
00941       fact=1.0;
00942       for (i=1;i<=MAXIT;i++)
00943       {
00944         fact *= -x/i;
00945         if (i != nm1)
00946         {
00947           del = -fact/(i-nm1);
00948         }
00949         else
00950         {
00951           psi = -EULER;  /* Compute psi(n) */
00952           for (ii=1;ii<=nm1;ii++)
00953           {
00954             psi += 1.0/ii;
00955           }
00956           del=fact*(-log(x)+psi);
00957         }
00958         ans += del;
00959         if (fabs(del) < fabs(ans)*EPS)
00960         {
00961           return ans;
00962         }
00963       }
00964       SevereError("expint", "series failed");
00965     }
00966   }
00967   return ans;
00968 }
00969 #undef EPS
00970 #undef EULER
00971 #undef MAXIT
00972 #undef FPMIN
00973 
00974 // ************ for adaptive quadrature ******************
00975 // ------------------- constructor ------------------------
00976 quad_list::quad_list( double ( *F )( double x, Param& params ) )
00977 {
00978   F_ = F;
00979 }
00980 // ------------------- quad_list::initialize ------------------
00981 // Set up the first link and get a rough approximate integral
00982 void quad_list::initialize( double A, double B, Param& params,
00983   double tol )
00984 {
00985   // start with a clean slate
00986   if( !empty( ) )
00987   {
00988     SevereError("quad_list::initialize", "Start with an empty list");
00989   }
00990 
00991   A_ = A;
00992   B_ = B;
00993   EPS_ = 2.0e-14;
00994   tol_ = ( tol < EPS_ ) ? EPS_ : tol;
00995   params_ = params;
00996 
00997   quad_link first_link;
00998 
00999   first_link.a = A_;
01000   first_link.b = B_;
01001   first_link.fa = F_( A_, params_ );
01002   first_link.fb = F_( B_, params_ );
01003   double x = 0.5*( A_ + B_ );
01004   first_link.fm = F_( x, params_ );
01005   insert( end( ), first_link );
01006 
01007   // for the approximate integral
01008   is_ = first_link.fa + first_link.fb + first_link.fm;
01009 
01010   // some "random" other points
01011   x = A_ + 0.9501 * ( B_ - A_ );
01012   is_ += F_( x, params_ );
01013 
01014   x = A_ + 0.2311 * ( B_ - A_ );
01015   is_ += F_( x, params_ );
01016 
01017   x = A_ + 0.6068 * ( B_ - A_ );
01018   is_ += F_( x, params_ );
01019 
01020   x = A_ + 0.4860 * ( B_ - A_ );
01021   is_ += F_( x, params_ );
01022 
01023   x = A_ + 0.8913 * ( B_ - A_ );
01024   is_ += F_( x, params_ );
01025 
01026   is_ *= ( B_ - A_ ) / 8;
01027 
01028   // ensure a nonzero approximation
01029   if( is_ == 0 )
01030   {
01031     is_ = B_ - A_;
01032   }
01033 
01034   // scale to account for the tolerance
01035   is_ *= tol_ / EPS_;
01036 
01037   // initialize the integral
01038   sum_ = 0.0;
01039 
01040   warning_set = false;
01041 }
01042 // ------------------- quad_list::test_int ------------------
01043 // Evaluate the integral over one interval and test its accuracy.
01044 void quad_list::test_int( quad_list::iterator link_ptr )
01045 {
01046   double m = ( link_ptr->a + link_ptr->b ) / 2;
01047   double h = ( link_ptr->b - link_ptr->a ) / 4;
01048   double fml = F_( link_ptr->a + h, params_ );
01049   double fmr = F_( link_ptr->b - h, params_ );
01050 
01051   // Simpson rule on original interval
01052   double i1 = ( link_ptr->fa + 4*link_ptr->fm + link_ptr->fb ) * h / 1.5;
01053 
01054   // Simpson rule on split interval
01055   double i2 = ( link_ptr->fa + 2*link_ptr->fm + link_ptr->fb +
01056                 4*( fml + fmr ) ) * h / 3;
01057 
01058   // Do a Richardson extrapolation
01059   i1 = ( 16*i2 - i1 ) / 15;
01060 
01061   // test the accuracy
01062   if( ( is_ + (i1 - i2 ) == is_ ) ||
01063       ( m <= link_ptr->a ) ||
01064       ( link_ptr->b <= m ) )
01065   {
01066     if( ( m <= link_ptr->a ) || ( link_ptr->b <= m ) )
01067     {
01068       warning_set = true;
01069     }
01070     sum_ += i1;  // update the running integral
01071     erase( link_ptr );  // we are done with this interval
01072   }
01073   else
01074   {
01075     // split this interval
01076     quad_link left_half;
01077     left_half.a = link_ptr->a;
01078     left_half.b = m;
01079     left_half.fa = link_ptr->fa;
01080     left_half.fb = link_ptr->fm;
01081     left_half.fm = fml;
01082     insert( link_ptr, left_half );
01083 
01084     // update the right half-interval
01085     link_ptr->a = m;
01086     link_ptr->fa = link_ptr->fm;
01087     link_ptr->fm = fmr;
01088   }
01089 }
01090 
01091 // ------------------- quad_list::Simp_quad ------------------
01092 // Evaluate the integral
01093 double quad_list::Simp_quad( double A, double B, Param& params,
01094   double tol )
01095 {
01096   int count = 0;
01097   initialize( A, B, params, tol );
01098 
01099   // process all of the subintervals
01100   while( !empty( ) )
01101   {
01102     quad_list::iterator link_ptr = begin( );
01103     // cout << "a: " << link_ptr->a << " b: " << link_ptr->b << endl;
01104     ++count;
01105     test_int( link_ptr );
01106   }
01107   if( warning_set )
01108   {
01109     Warning("quad_list::Simp_quad", "The integral may not be accurate");
01110   }
01111   // cout << "intervals: " << count << endl;
01112   return sum_;
01113 }

Generated on Thu Sep 7 10:30:07 2006 for fete -- From ENDFB6 To ENDL by doxygen 1.3.4