# include # include # include # include # include # include # include # include using namespace std; int main ( int argc, char *argv[] ); char ch_cap ( char ch ); bool ch_eqi ( char ch1, char ch2 ); int ch_to_digit ( char ch ); void dtable_close_write ( ofstream &output ); void dtable_data_write ( ofstream &output, int m, int n, double table[] ); void dtable_header_write ( char *output_filename, ofstream &output, int m, int n ); void dtable_write ( char *output_filename, int m, int n, double table[], bool header ); void gegenbauer_compute ( int order, double alpha, double xtab[], double weight[] ); void gegenbauer_handle ( int order, double alpha, char *output ); void gegenbauer_recur ( double *p2, double *dp2, double *p1, double x, int order, double alpha, double b[], double c[] ); void gegenbauer_root ( double *x, int order, double alpha, double *dp2, double *p1, double c[] ); double r8_abs ( double x ); double r8_epsilon ( void ); double r8_gamma ( double x ); double r8_huge ( void ); bool s_eqi ( char *s1, char *s2 ); void timestamp ( void ); char *timestring ( void ); //****************************************************************************80 int main ( int argc, char *argv[] ) //****************************************************************************80 // // Purpose: // // MAIN is the main program for GEGENBAUER_RULE. // // Discussion: // // This program computes a standard Gauss-Gegenbauer quadrature rule // and writes it to a file. // // The user specifies: // * the ORDER (number of points) in the rule // * the ALPHA parameter, // * the root name of the output files. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 01 March 2008 // // Author: // // John Burkardt // { double alpha; int order; char output[81]; cout << "\n"; timestamp ( ); cout << "\n"; cout << "GEGENBAUER_RULE\n"; cout << " C++ version\n"; cout << "\n"; cout << " Compiled on " << __DATE__ << " at " << __TIME__ << ".\n"; cout << "\n"; cout << " Compute a Gauss-Gegenbauer quadrature rule for approximating\n"; cout << "\n"; cout << " Integral ( -1 <= x <= +1 ) (1-x^2)^ALPHA f(x) dx\n"; cout << "\n"; cout << " of order ORDER.\n"; cout << "\n"; cout << " The user specifies ORDER, ALPHA, and OUTPUT.\n"; cout << "\n"; cout << " OUTPUT is:\n"; cout << "\n"; cout << " \"C++\" for printed C++ output;\n"; cout << " \"F77\" for printed Fortran77 output;\n"; cout << " \"F90\" for printed Fortran90 output;\n"; cout << " \"MAT\" for printed MATLAB output;\n"; cout << "\n"; cout << " or:\n"; cout << "\n"; cout << " \"filename\" to generate 3 files:\n"; cout << "\n"; cout << " filename_w.txt - the weight file\n"; cout << " filename_x.txt - the abscissa file.\n"; cout << " filename_r.txt - the region file.\n"; // // Get the order. // if ( 1 < argc ) { order = atoi ( argv[1] ); } else { cout << "\n"; cout << " Enter the value of ORDER (1 or greater)\n"; cin >> order; } cout << "\n"; cout << " The requested order of the rule is = " << order << "\n"; // // Get ALPHA. // if ( 2 < argc ) { alpha = atof ( argv[2] ); } else { cout << "\n"; cout << " ALPHA is the exponent of (1-x^2) in the integral:\n"; cout << "\n"; cout << " Integral ( -1 <= x <= +1 ) (1-x^2)^ALPHA f(x) dx\n"; cout << "\n"; cout << " Note that -1.0 < ALPHA is required.\n"; cout << "\n"; cout << " Enter the value of ALPHA:\n"; cin >> alpha; } cout << "\n"; cout << " The requested value of ALPHA = " << alpha << "\n"; // // Get the output option or quadrature file root name: // if ( 3 < argc ) { strcpy ( output, argv[3] ); } else { cout << "\n"; cout << " Enter OUTPUT (one of C++, F77, F90, MAT\n"; cout << " or else the \"root name\" of the quadrature files).\n"; cin >> output; } cout << "\n"; cout << " OUTPUT option is \"" << output << "\".\n"; // // Construct the rule and output it. // gegenbauer_handle ( order, alpha, output ); cout << "\n"; cout << "GEGENBAUER_RULE:\n"; cout << " Normal end of execution.\n"; cout << "\n"; timestamp ( ); return 0; } //****************************************************************************80 char ch_cap ( char ch ) //****************************************************************************80 // // Purpose: // // CH_CAP capitalizes a single character. // // Discussion: // // This routine should be equivalent to the library "toupper" function. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 19 July 1998 // // Author: // // John Burkardt // // Parameters: // // Input, char CH, the character to capitalize. // // Output, char CH_CAP, the capitalized character. // { if ( 97 <= ch && ch <= 122 ) { ch = ch - 32; } return ch; } //****************************************************************************80 bool ch_eqi ( char ch1, char ch2 ) //****************************************************************************80 // // Purpose: // // CH_EQI is true if two characters are equal, disregarding case. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 13 June 2003 // // Author: // // John Burkardt // // Parameters: // // Input, char CH1, CH2, the characters to compare. // // Output, bool CH_EQI, is true if the two characters are equal, // disregarding case. // { if ( 97 <= ch1 && ch1 <= 122 ) { ch1 = ch1 - 32; } if ( 97 <= ch2 && ch2 <= 122 ) { ch2 = ch2 - 32; } return ( ch1 == ch2 ); } //****************************************************************************80 int ch_to_digit ( char ch ) //****************************************************************************80 // // Purpose: // // CH_TO_DIGIT returns the integer value of a base 10 digit. // // Example: // // CH DIGIT // --- ----- // '0' 0 // '1' 1 // ... ... // '9' 9 // ' ' 0 // 'X' -1 // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 13 June 2003 // // Author: // // John Burkardt // // Parameters: // // Input, char CH, the decimal digit, '0' through '9' or blank are legal. // // Output, int CH_TO_DIGIT, the corresponding integer value. If the // character was 'illegal', then DIGIT is -1. // { int digit; if ( '0' <= ch && ch <= '9' ) { digit = ch - '0'; } else if ( ch == ' ' ) { digit = 0; } else { digit = -1; } return digit; } //****************************************************************************80 void dtable_close_write ( ofstream &output ) //****************************************************************************80 // // Purpose: // // DTABLE_CLOSE_WRITE closes a file to which a DTABLE was to be written. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 12 November 2007 // // Author: // // John Burkardt // // Parameters: // // Input, ofstream output, the output file stream. // { output.close ( ); return; } //****************************************************************************80 void dtable_data_write ( ofstream &output, int m, int n, double table[] ) //****************************************************************************80 // // Purpose: // // DTABLE_DATA_WRITE writes data to a DTABLE file. // // Discussion: // // The file should already be open. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 12 November 2007 // // Author: // // John Burkardt // // Parameters: // // Input, ofstream &OUTPUT, a pointer to the output stream. // // Input, int M, the spatial dimension. // // Input, int N, the number of points. // // Input, double TABLE[M*N], the table data. // { int i; int j; char *s; for ( j = 0; j < n; j++ ) { for ( i = 0; i < m; i++ ) { output << setprecision(16) << setw(24) << table[i+j*m] << " "; } output << "\n"; } output.close ( ); return; } //****************************************************************************80 void dtable_header_write ( char *output_filename, ofstream &output, int m, int n ) //****************************************************************************80 // // Purpose: // // DTABLE_HEADER_WRITE writes the header of a DTABLE file. // // Discussion: // // The file should already be open. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 12 November 2007 // // Author: // // John Burkardt // // Parameters: // // Input, char *OUTPUT_FILENAME, the output filename. // // Input, ofstream &OUTPUT, the output stream. // // Input, int M, the spatial dimension. // // Input, int N, the number of points. // { char *s; s = timestring ( ); output << "# " << output_filename << "\n"; output << "# created by dtable_write in table_io.C" << "\n"; output << "# at " << s << "\n"; output << "#\n"; output << "# Spatial dimension M = " << m << "\n"; output << "# Number of points N = " << n << "\n"; output << "# EPSILON (unit roundoff) = " << r8_epsilon ( ) << "\n"; output << "#\n"; delete [] s; return; } //****************************************************************************80 void dtable_write ( char *output_filename, int m, int n, double table[], bool header ) //****************************************************************************80 // // Purpose: // // DTABLE_WRITE writes information to a DTABLE file. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 12 November 2007 // // Author: // // John Burkardt // // Parameters: // // Input, char *OUTPUT_FILENAME, the output filename. // // Input, int M, the spatial dimension. // // Input, int N, the number of points. // // Input, double TABLE[M*N], the table data. // // Input, bool HEADER, is TRUE if the header is to be included. // { ofstream output; output.open ( output_filename ); if ( !output ) if ( false ) { cout << "\n"; cout << "DTABLE_WRITE - Fatal error!\n"; cout << " Could not open the output file.\n"; return; } if ( header ) { dtable_header_write ( output_filename, output, m, n ); } dtable_data_write ( output, m, n, table ); dtable_close_write ( output ); return; } //****************************************************************************80 void gegenbauer_compute ( int order, double alpha, double xtab[], double weight[] ) //****************************************************************************80 // // Purpose: // // GEGENBAUER_COMPUTE computes a Gauss-Gegenbauer quadrature rule. // // Discussion: // // The integration interval is [ -1, 1 ]. // // The integral to approximate: // // Integral ( -1 <= X <= 1 ) (1-X^2)^ALPHA * F(X) dX // // The quadrature rule: // // Sum ( 1 <= I <= ORDER ) WEIGHT(I) * F ( XTAB(I) ) // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 26 February 2008 // // Author: // // John Burkardt // // Reference: // // Arthur Stroud, Don Secrest, // Gaussian Quadrature Formulas, // Prentice Hall, 1966, // LC: QA299.4G3S7. // // Parameters: // // Input, int ORDER, the order of the quadrature rule to be computed. // // Input, double ALPHA, the exponent of (1-X^2). -1.0 < ALPHA is required. // // Output, double XTAB[ORDER], the abscissas. // // Output, double WEIGHT[ORDER], the weights. // { double an; double *c; double cc; double delta; double dp2; int i; double p1; double prod; double r1; double r2; double r3; double temp; double x; c = new double[order]; // // Check ALPHA. // if ( alpha <= -1.0 ) { cout << "\n"; cout << "GEGENBAUER_COMPUTE - Fatal error!\n"; cout << " -1.0 < ALPHA is required.\n"; exit ( 1 ); } // // Set the recursion coefficients. // for ( i = 1; i <= order; i++ ) { if ( i == 1 ) { c[i-1] = 0.0; } else { c[i-1] = ( double ) ( i - 1 ) * ( alpha + alpha + ( double ) ( i - 1 ) ) / ( ( alpha + alpha + ( double ) ( 2 * i - 1 ) ) * ( alpha + alpha + ( double ) ( 2 * i - 3 ) ) ); } } delta = r8_gamma ( alpha + 1.0 ) * r8_gamma ( alpha + 1.0 ) / r8_gamma ( alpha + alpha + 2.0 ); prod = 1.0; for ( i = 2; i <= order; i++ ) { prod = prod * c[i-1]; } cc = delta * pow ( 2.0, alpha + alpha + 1.0 ) * prod; for ( i = 1; i <= order; i++ ) { if ( i == 1 ) { an = alpha / ( double ) ( order ); r1 = ( 1.0 + alpha ) * ( 2.78 / ( 4.0 + ( double ) ( order * order ) ) + 0.768 * an / ( double ) ( order ) ); r2 = 1.0 + 2.44 * an + 1.282 * an * an; x = ( r2 - r1 ) / r2; } else if ( i == 2 ) { r1 = ( 4.1 + alpha ) / ( ( 1.0 + alpha ) * ( 1.0 + 0.156 * alpha ) ); r2 = 1.0 + 0.06 * ( ( double ) ( order ) - 8.0 ) * ( 1.0 + 0.12 * alpha ) / ( double ) ( order ); r3 = 1.0 + 0.012 * alpha * ( 1.0 + 0.25 * r8_abs ( alpha ) ) / ( double ) ( order ); x = x - r1 * r2 * r3 * ( 1.0 - x ); } else if ( i == 3 ) { r1 = ( 1.67 + 0.28 * alpha ) / ( 1.0 + 0.37 * alpha ); r2 = 1.0 + 0.22 * ( ( double ) ( order ) - 8.0 ) / ( double ) ( order ); r3 = 1.0 + 8.0 * alpha / ( ( 6.28 + alpha ) * ( double ) ( order * order ) ); x = x - r1 * r2 * r3 * ( xtab[0] - x ); } else if ( i < order - 1 ) { x = 3.0 * xtab[i-2] - 3.0 * xtab[i-3] + xtab[i-4]; } else if ( i == order - 1 ) { r1 = ( 1.0 + 0.235 * alpha ) / ( 0.766 + 0.119 * alpha ); r2 = 1.0 / ( 1.0 + 0.639 * ( ( double ) ( order ) - 4.0 ) / ( 1.0 + 0.71 * ( ( double ) ( order ) - 4.0 ) ) ); r3 = 1.0 / ( 1.0 + 20.0 * alpha / ( ( 7.5 + alpha ) * ( double ) ( order * order ) ) ); x = x + r1 * r2 * r3 * ( x - xtab[i-3] ); } else if ( i == order ) { r1 = ( 1.0 + 0.37 * alpha ) / ( 1.67 + 0.28 * alpha ); r2 = 1.0 / ( 1.0 + 0.22 * ( ( double ) ( order ) - 8.0 ) / ( double ) ( order ) ); r3 = 1.0 / ( 1.0 + 8.0 * alpha / ( ( 6.28 + alpha ) * ( double ) ( order * order ) ) ); x = x + r1 * r2 * r3 * ( x - xtab[i-3] ); } gegenbauer_root ( &x, order, alpha, &dp2, &p1, c ); xtab[i-1] = x; weight[i-1] = cc / ( dp2 * p1 ); } // // Reverse the order of the values. // for ( i = 1; i <= order/2; i++ ) { temp = xtab[i-1]; xtab[i-1] = xtab[order-i]; xtab[order-i] = temp; } for ( i = 1; i <=order/2; i++ ) { temp = weight[i-1]; weight[i-1] = weight[order-i]; weight[order-i] = temp; } delete [] c; return; } //****************************************************************************80 void gegenbauer_handle ( int order, double alpha, char *output ) //****************************************************************************80 // // Purpose: // // GEGENBAUER_HANDLE handles the Gauss-Gegenbauer rule. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 01 March 2008 // // Author: // // John Burkardt // // Parameters: // // Input, int ORDER, the order of the rule. // // Input, double ALPHA, the parameter. // // Input, char *OUTPUT_FILE, specifies the output. // * "C++'", print as C++ code. // * "F77", print as FORTRAN77 code. // * "F90", print as FORTRAN90 code. // * "MAT", print as MATLAB code. // * file, write files "file_w.txt", "file_x.txt", "file_r.txt" // defining weights, abscissas, and region. // { bool header; int i; char output_r[81]; char output_w[81]; char output_x[81]; double *r; double *w; double *x; r = new double[2]; w = new double[order]; x = new double[order]; r[0] = - 1.0; r[1] = + 1.0; gegenbauer_compute ( order, alpha, x, w ); if ( s_eqi ( output, "C++" ) ) { cout << "//\n"; cout << "// Weights W, abscissas X and range R\n"; cout << "// for a Gauss-Gegenbauer quadrature rule\n"; cout << "// ORDER = " << order << "\n"; cout << "// ALPHA = " << alpha << "\n"; cout << "//\n"; cout << "// Standard rule:\n"; cout << "// Integral ( -1 <= x <= +1 ) (1-x^2)^ALPHA f(x) dx\n"; cout << "// is to be approximated by\n"; cout << "// sum ( 1 <= I <= ORDER ) w(i) * f(x(i)).\n"; cout << "//\n"; for ( i = 0; i < order; i++ ) { cout << " w[" << i << "] = " << setprecision(16) << w[i] << ";\n"; } cout << "\n"; for ( i = 0; i < order; i++ ) { cout << " x[" << i << "] = " << setprecision(16) << x[i] << ";\n"; } cout << "\n"; for ( i = 0; i < 2; i++ ) { cout << " r[" << i << "] = " << r[i] << ";\n"; } } else if ( s_eqi ( output, "F77" ) ) { cout << "c\n"; cout << "c Weights W, abscissas X and range R\n"; cout << "c for a Gauss-Gegenbauer quadrature rule\n"; cout << "c ORDER = " << order << "\n"; cout << "c ALPHA = " << alpha << "\n"; cout << "c\n"; cout << "c Standard rule:\n"; cout << "c Integral ( -1 <= x <= +1 ) (1-x^2)^ALPHA f(x) dx\n"; cout << "c is to be approximated by\n"; cout << "c sum ( 1 <= I <= ORDER ) w(i) * f(x(i)).\n"; cout << "c\n"; for ( i = 0; i < order; i++ ) { cout << " w(" << i + 1 << ") = " << setprecision(16) << w[i] << "\n"; } cout << "\n"; for ( i = 0; i < order; i++ ) { cout << " x(" << i + 1 << ") = " << setprecision(16) << x[i] << "\n"; } cout << "\n"; for ( i = 0; i < 2; i++ ) { cout << " r(" << i + 1 << ") = " << r[i] << "\n"; } } else if ( s_eqi ( output, "F90" ) ) { cout << "!\n"; cout << "! Weights W, abscissas X and range R\n"; cout << "! for a Gauss-Gegenbauer quadrature rule\n"; cout << "! ORDER = " << order << "\n"; cout << "! ALPHA = " << alpha << "\n"; cout << "!\n"; cout << "! Standard rule:\n"; cout << "! Integral ( -1 <= x <= +1 ) (1-x^2)^ALPHA f(x) dx\n"; cout << "! is to be approximated by\n"; cout << "! sum ( 1 <= I <= ORDER ) w(i) * f(x(i)).\n"; cout << "!\n"; for ( i = 0; i < order; i++ ) { cout << " w(" << i + 1 << ") = " << setprecision(16) << w[i] << "\n"; } cout << "\n"; for ( i = 0; i < order; i++ ) { cout << " x(" << i + 1 << ") = " << setprecision(16) << x[i] << "\n"; } cout << "\n"; for ( i = 0; i < 2; i++ ) { cout << " r(" << i + 1 << ") = " << r[i] << "\n"; } } else if ( s_eqi ( output, "MAT" ) ) { cout << "%\n"; cout << "% Weights W, abscissas X and range R\n"; cout << "% for a Gauss-Gegenbauer quadrature rule\n"; cout << "% ORDER = " << order << "\n"; cout << "% ALPHA = " << alpha << "\n"; cout << "%\n"; cout << "% Standard rule:\n"; cout << "% Integral ( -1 <= x <= +1 ) (1-x^2)^ALPHA f(x) dx\n"; cout << "% is to be approximated by\n"; cout << "% sum ( 1 <= I <= ORDER ) w(i) * f(x(i)).\n"; cout << "%\n"; for ( i = 0; i < order; i++ ) { cout << " w(" << i + 1 << ") = " << setprecision(16) << w[i] << ";\n"; } cout << "\n"; for ( i = 0; i < order; i++ ) { cout << " x(" << i + 1 << ") = " << setprecision(16) << x[i] << ";\n"; } cout << "\n"; for ( i = 0; i < 2; i++ ) { cout << " r(" << i + 1 << ") = " << r[i] << ";\n"; } } else { sprintf ( output_w, "%s_w.txt", output ); sprintf ( output_x, "%s_x.txt", output ); sprintf ( output_r, "%s_r.txt", output ); cout << "\n"; cout << " Creating quadrature files.\n"; cout << "\n"; cout << " Root file name is \"" << output << "\".\n"; cout << "\n"; cout << " Weight file will be \"" << output_w << "\".\n"; cout << " Abscissa file will be \"" << output_x << "\".\n"; cout << " Region file will be \"" << output_r << "\".\n"; header = false; dtable_write ( output_w, 1, order, w, header ); dtable_write ( output_x, 1, order, x, header ); dtable_write ( output_r, 1, 2, r, header ); } delete [] r; delete [] w; delete [] x; return; } //****************************************************************************80 void gegenbauer_recur ( double *p2, double *dp2, double *p1, double x, int order, double alpha, double c[] ) //****************************************************************************80 // // Purpose: // // GEGENBAUER_RECUR finds the value and derivative of a Gegenbauer polynomial. // // Modified: // // 26 February 2008 // // Author: // // John Burkardt // // Reference: // // Arthur Stroud, Don Secrest, // Gaussian Quadrature Formulas, // Prentice Hall, 1966, // LC: QA299.4G3S7. // // Parameters: // // Output, double *P2, the value of J(ORDER)(X). // // Output, double *DP2, the value of J'(ORDER)(X). // // Output, double *P1, the value of J(ORDER-1)(X). // // Input, double X, the point at which polynomials are evaluated. // // Input, int ORDER, the order of the polynomial to be computed. // // Input, double ALPHA, the exponents of (1-X^2). // // Input, double C[ORDER], the recursion coefficients. // { double dp0; double dp1; int i; double p0; *p1 = 1.0; dp1 = 0.0; *p2 = x; *dp2 = 1.0; for ( i = 2; i <= order; i++ ) { p0 = *p1; dp0 = dp1; *p1 = *p2; dp1 = *dp2; *p2 = x * ( *p1 ) - c[i-1] * p0; *dp2 = x * dp1 + ( *p1 ) - c[i-1] * dp0; } return; } //****************************************************************************80 void gegenbauer_root ( double *x, int order, double alpha, double *dp2, double *p1, double c[] ) //****************************************************************************80 // // Purpose: // // GEGENBAUER_ROOT improves an approximate root of a Gegenbauer polynomial. // // Modified: // // 26 February 2008 // // Author: // // John Burkardt // // Reference: // // Arthur Stroud, Don Secrest, // Gaussian Quadrature Formulas, // Prentice Hall, 1966, // LC: QA299.4G3S7. // // Parameters: // // Input/output, double *X, the approximate root, which // should be improved on output. // // Input, int ORDER, the order of the polynomial to be computed. // // Input, double ALPHA, the exponents of (1-X^2). // // Output, double *DP2, the value of J'(ORDER)(X). // // Output, double *P1, the value of J(ORDER-1)(X). // // Input, double C[ORDER], the recursion coefficients. // { double d; double eps; double p2; int step; int step_max = 10; eps = r8_epsilon ( ); for ( step = 1; step <= step_max; step++ ) { gegenbauer_recur ( &p2, dp2, p1, *x, order, alpha, c ); d = p2 / ( *dp2 ); *x = *x - d; if ( r8_abs ( d ) <= eps * ( r8_abs ( *x ) + 1.0 ) ) { return; } } return; } //****************************************************************************80 double r8_abs ( double x ) //****************************************************************************80 // // Purpose: // // R8_ABS returns the absolute value of an R8. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 14 November 2006 // // Author: // // John Burkardt // // Parameters: // // Input, double X, the quantity whose absolute value is desired. // // Output, double R8_ABS, the absolute value of X. // { double value; if ( 0.0 <= x ) { value = x; } else { value = -x; } return value; } //****************************************************************************80 double r8_epsilon ( void ) //****************************************************************************80 // // Purpose: // // R8_EPSILON returns the R8 roundoff unit. // // Discussion: // // The roundoff unit is a number R which is a power of 2 with the // property that, to the precision of the computer's arithmetic, // 1 < 1 + R // but // 1 = ( 1 + R / 2 ) // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 01 July 2004 // // Author: // // John Burkardt // // Parameters: // // Output, double R8_EPSILON, the R8 round-off unit. // { double value; value = 1.0; while ( 1.0 < ( double ) ( 1.0 + value ) ) { value = value / 2.0; } value = 2.0 * value; return value; } //****************************************************************************80 double r8_gamma ( double x ) //****************************************************************************80 // // Purpose: // // R8_GAMMA evaluates Gamma(X) for a real argument. // // Discussion: // // This function was originally named DGAMMA. // // However, a number of Fortran compilers now include a library // function of this name. To avoid conflicts, this function was // renamed R8_GAMMA. // // This routine calculates the GAMMA function for a real argument X. // Computation is based on an algorithm outlined in reference 1. // The program uses rational functions that approximate the GAMMA // function to at least 20 significant decimal digits. Coefficients // for the approximation over the interval (1,2) are unpublished. // Those for the approximation for 12 <= X are from reference 2. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 18 January 2008 // // Author: // // FORTRAN77 original version by William Cody, Laura Stoltz // C++ version by John Burkardt // // Reference: // // William Cody, // An Overview of Software Development for Special Functions, // in Numerical Analysis Dundee, 1975, // edited by GA Watson, // Lecture Notes in Mathematics 506, // Springer, 1976. // // John Hart, Ward Cheney, Charles Lawson, Hans Maehly, // Charles Mesztenyi, John Rice, Henry Thatcher, // Christoph Witzgall, // Computer Approximations, // Wiley, 1968, // LC: QA297.C64. // // Parameters: // // Input, double X, the argument of the function. // // Output, double R8_GAMMA, the value of the function. // { // // Coefficients for minimax approximation over (12, INF). // double c[7] = { -1.910444077728E-03, 8.4171387781295E-04, -5.952379913043012E-04, 7.93650793500350248E-04, -2.777777777777681622553E-03, 8.333333333333333331554247E-02, 5.7083835261E-03 }; double eps = 2.22E-16; double fact; double half = 0.5; int i; int n; double one = 1.0; double p[8] = { -1.71618513886549492533811E+00, 2.47656508055759199108314E+01, -3.79804256470945635097577E+02, 6.29331155312818442661052E+02, 8.66966202790413211295064E+02, -3.14512729688483675254357E+04, -3.61444134186911729807069E+04, 6.64561438202405440627855E+04 }; bool parity; double pi = 3.1415926535897932384626434; double q[8] = { -3.08402300119738975254353E+01, 3.15350626979604161529144E+02, -1.01515636749021914166146E+03, -3.10777167157231109440444E+03, 2.25381184209801510330112E+04, 4.75584627752788110767815E+03, -1.34659959864969306392456E+05, -1.15132259675553483497211E+05 }; double res; double sqrtpi = 0.9189385332046727417803297; double sum; double twelve = 12.0; double two = 2.0; double value; double xbig = 171.624; double xden; double xinf = 1.79E+308; double xminin = 2.23E-308; double xnum; double y; double y1; double ysq; double z; double zero = 0.0;; parity = false; fact = one; n = 0; y = x; // // Argument is negative. // if ( y <= zero ) { y = - x; y1 = ( double ) ( int ) ( y ); res = y - y1; if ( res != zero ) { if ( y1 != ( double ) ( int ) ( y1 * half ) * two ) { parity = true; } fact = - pi / sin ( pi * res ); y = y + one; } else { res = xinf; value = res; return value; } } // // Argument is positive. // if ( y < eps ) { // // Argument < EPS. // if ( xminin <= y ) { res = one / y; } else { res = xinf; value = res; return value; } } else if ( y < twelve ) { y1 = y; // // 0.0 < argument < 1.0. // if ( y < one ) { z = y; y = y + one; } // // 1.0 < argument < 12.0. // Reduce argument if necessary. // else { n = ( int ) ( y ) - 1; y = y - ( double ) ( n ); z = y - one; } // // Evaluate approximation for 1.0 < argument < 2.0. // xnum = zero; xden = one; for ( i = 0; i < 8; i++ ) { xnum = ( xnum + p[i] ) * z; xden = xden * z + q[i]; } res = xnum / xden + one; // // Adjust result for case 0.0 < argument < 1.0. // if ( y1 < y ) { res = res / y1; } // // Adjust result for case 2.0 < argument < 12.0. // else if ( y < y1 ) { for ( i = 1; i <= n; i++ ) { res = res * y; y = y + one; } } } else { // // Evaluate for 12.0 <= argument. // if ( y <= xbig ) { ysq = y * y; sum = c[6]; for ( i = 0; i < 6; i++ ) { sum = sum / ysq + c[i]; } sum = sum / y - y + sqrtpi; sum = sum + ( y - half ) * log ( y ); res = exp ( sum ); } else { res = xinf; value = res; return value; } } // // Final adjustments and return. // if ( parity ) { res = - res; } if ( fact != one ) { res = fact / res; } value = res; return value; } //****************************************************************************80 double r8_huge ( void ) //****************************************************************************80 // // Purpose: // // R8_HUGE returns a "huge" R8. // // Discussion: // // The value returned by this function is NOT required to be the // maximum representable R8. This value varies from machine to machine, // from compiler to compiler, and may cause problems when being printed. // We simply want a "very large" but non-infinite number. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 06 October 2007 // // Author: // // John Burkardt // // Parameters: // // Output, double R8_HUGE, a "huge" R8 value. // { double value; value = 1.0E+30; return value; } //****************************************************************************80 bool s_eqi ( char *s1, char *s2 ) //****************************************************************************80 // // Purpose: // // S_EQI reports whether two strings are equal, ignoring case. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 05 May 2003 // // Author: // // John Burkardt // // Parameters: // // Input, char *S1, char *S2, pointers to two strings. // // Output, bool S_EQI, is true if the strings are equal. // { int i; int nchar; int nchar1; int nchar2; nchar1 = strlen ( s1 ); nchar2 = strlen ( s2 ); if ( nchar1 < nchar2 ) { nchar = nchar1; } else { nchar = nchar2; } // // The strings are not equal if they differ over their common length. // for ( i = 0; i < nchar; i++ ) { if ( ch_cap ( s1[i] ) != ch_cap ( s2[i] ) ) { return false; } } // // The strings are not equal if the longer one includes nonblanks // in the tail. // if ( nchar < nchar1 ) { for ( i = nchar; i < nchar1; i++ ) { if ( s1[i] != ' ' ) { return false; } } } else if ( nchar < nchar2 ) { for ( i = nchar; i < nchar2; i++ ) { if ( s2[i] != ' ' ) { return false; } } } return true; } //****************************************************************************80 void timestamp ( void ) //****************************************************************************80 // // Purpose: // // TIMESTAMP prints the current YMDHMS date as a time stamp. // // Example: // // 31 May 2001 09:45:54 AM // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 24 September 2003 // // Author: // // John Burkardt // // Parameters: // // None // { # define TIME_SIZE 40 static char time_buffer[TIME_SIZE]; const struct tm *tm; size_t len; time_t now; now = time ( NULL ); tm = localtime ( &now ); len = strftime ( time_buffer, TIME_SIZE, "%d %B %Y %I:%M:%S %p", tm ); cout << time_buffer << "\n"; return; # undef TIME_SIZE } //****************************************************************************80 char *timestring ( void ) //****************************************************************************80 // // Purpose: // // TIMESTRING returns the current YMDHMS date as a string. // // Example: // // 31 May 2001 09:45:54 AM // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 24 September 2003 // // Author: // // John Burkardt // // Parameters: // // Output, char *TIMESTRING, a string containing the current YMDHMS date. // { # define TIME_SIZE 40 const struct tm *tm; size_t len; time_t now; char *s; now = time ( NULL ); tm = localtime ( &now ); s = new char[TIME_SIZE]; len = strftime ( s, TIME_SIZE, "%d %B %Y %I:%M:%S %p", tm ); return s; # undef TIME_SIZE }