# include # include # include # include # include # include using namespace std; # include "gegenbauer_polynomial.hpp" int main ( ); void gegenbauer_alpha_check_test ( ); void gegenbauer_ek_compute_test ( ); void gegenbauer_integral_test ( ); void gegenbauer_polynomial_value_test ( ); void gegenbauer_ss_compute_test ( ); void gegenbauer_to_monomial_matrix_test ( ); void imtqlx_test ( ); void monomial_to_gegenbauer_matrix_test ( ); void r8_hyper_2f1_test ( ); void r8_uniform_ab_test ( ); double *r8mat_mm_new ( int n1, int n2, int n3, double a[], double b[] ); double *r8mat_mv_new ( int m, int n, double a[], double x[] ); void r8mat_print ( int m, int n, double a[], string title ); void r8mat_print_some ( int m, int n, double a[], int ilo, int jlo, int ihi, int jhi, string title ); void r8poly_print ( int n, double a[], string title ); void timestamp ( ); //****************************************************************************80 int main ( ) //****************************************************************************80 // // Purpose: // // gegenbauer_polynomial_test() tests gegenbauer_polynomial(). // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 24 April 2024 // // Author: // // John Burkardt // { timestamp ( ); cout << "\n"; cout << "gegenbauer_polynomial_test():\n"; cout << " C++ version.\n"; cout << " Test gegenbauer_polynomial().\n"; gegenbauer_alpha_check_test ( ); gegenbauer_ek_compute_test ( ); gegenbauer_integral_test ( ); gegenbauer_polynomial_value_test ( ); gegenbauer_ss_compute_test ( ); gegenbauer_to_monomial_matrix_test ( ); imtqlx_test ( ); monomial_to_gegenbauer_matrix_test ( ); r8_hyper_2f1_test ( ); r8_uniform_ab_test ( ); // // Terminate. // cout << "\n"; cout << "gegenbauer_polynomial_test():\n"; cout << " Normal end of execution.\n"; cout << "\n"; timestamp ( ); return 0; } //****************************************************************************80 void gegenbauer_alpha_check_test ( ) //****************************************************************************80 // // Purpose: // // gegenbauer_alpha_check_test() tests gegenbauer_alpha_check(). // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 30 November 2015 // // Author: // // John Burkardt // { double alpha; bool check; int n; int seed; cout << "\n"; cout << "gegenbauer_alpha_check_test():\n"; cout << " gegenbauer_alpha_check() checks that ALPHA is legal;\n"; cout << "\n"; cout << " ALPHA Check?\n"; cout << "\n"; seed = 123456789; for ( n = 1; n <= 10; n++ ) { alpha = r8_uniform_ab ( -5.0, +5.0, seed ); check = gegenbauer_alpha_check ( alpha ); cout << " " << setw(10) << alpha << " " << setw(1) << check << "\n"; } return; } //****************************************************************************80 void gegenbauer_ek_compute_test ( ) //****************************************************************************80 // // Purpose: // // gegenbauer_ek_compute_test() tests gegenbauer_ek_compute(). // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 20 November 2015 // // Author: // // John Burkardt // { double alpha; int i; int n; int prec; double *w; double *x; prec = cout.precision ( ); cout.precision ( 16 ); alpha = 0.5; cout << "\n"; cout << "gegenbauer_ek_compute_test():\n"; cout << " gegenbauer_ek_compute() computes a Gauss-Gegenbauer rule;\n"; cout << "\n"; cout << " with ALPHA = " << alpha << "\n"; cout << " and integration interval [-1,+1]\n"; cout << "\n"; cout << " W X\n"; for ( n = 1; n <= 10; n++ ) { cout << "\n"; w = new double[n]; x = new double[n]; gegenbauer_ek_compute ( n, alpha, x, w ); for ( i = 0; i < n; i++ ) { cout << " " << " " << setw(14) << w[i] << " " << setw(14) << x[i] << "\n"; } delete [] w; delete [] x; } cout.precision ( prec ); return; } //****************************************************************************80 void gegenbauer_integral_test ( ) //****************************************************************************80 // // Purpose: // // gegenbauer_integral_test() tests gegenbauer_integral(). // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 14 June 2015 // // Author: // // John Burkardt // { double alpha; int n; int prec; double value; prec = cout.precision ( ); cout.precision ( 16 ); alpha = 0.25; cout << "\n"; cout << "GEGENBAUER_INTEGRAL_TEST\n"; cout << " GEGENBAUER_INTEGRAL evaluates\n"; cout << " Integral ( -1 < x < +1 ) x^n * (1-x*x)^alpha dx\n"; cout << "\n"; cout << " N Value\n"; cout << "\n"; for ( n = 0; n <= 10; n++ ) { value = gegenbauer_integral ( n, alpha ); cout << " " << setw(8) << n << " " << setw(24) << value << "\n"; } cout.precision ( prec ); return; } //****************************************************************************80 void gegenbauer_polynomial_value_test ( ) //****************************************************************************80 // // Purpose: // // GEGENBAUER_POLYNOMIAL_VALUE_TEST tests GEGENBAUER_POLYNOMIAL_VALUE. // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 30 November 2015 // // Author: // // John Burkardt // { double alpha; double *c; double fx; int m; int n; int n_data; double x[1]; double xscalar; cout << "\n"; cout << "GEGENBAUER_POLYNOMIAL_VALUE_TEST:\n"; cout << " GEGENBAUER_POLYNOMIAL_VALUE evaluates the Gegenbauer polynomial.\n"; cout << "\n"; cout << " M ALPHA X GPV GEGENBAUER\n"; cout << "\n"; n_data = 0; for ( ; ; ) { gegenbauer_polynomial_values ( n_data, m, alpha, xscalar, fx ); if ( n_data == 0 ) { break; } // // Since GEGENBAUER_POLYNOMIAL_VALUE expects a vector input X, we have to // do a little "rewrapping" of the input. // n = 1; x[0] = xscalar; c = gegenbauer_polynomial_value ( m, n, alpha, x ); cout << " " << setw(6) << m << " " << setw(8) << alpha << " " << setw(8) << x[0] << " " << setw(12) << fx << " " << setw(12) << c[m+0*(m+1)] << "\n"; delete [] c; } return; } //****************************************************************************80 void gegenbauer_ss_compute_test ( ) //****************************************************************************80 // // Purpose: // // gegenbauer_ss_compute_test() tests gegenbauer_ss_compute(). // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 10 June 2015 // // Author: // // John Burkardt // { double alpha; int i; int n; int prec; double *w; double *x; prec = cout.precision ( ); cout.precision ( 16 ); alpha = 0.5; cout << "\n"; cout << "gegenbauer_ss_compute_test():\n"; cout << " gegenbauer_ss_compute() computes a Gauss-Gegenbauer rule;\n"; cout << "\n"; cout << " with ALPHA = " << alpha << "\n"; cout << "\n"; cout << " W X\n"; for ( n = 1; n <= 10; n++ ) { cout << "\n"; w = new double[n]; x = new double[n]; gegenbauer_ss_compute ( n, alpha, x, w ); for ( i = 0; i < n; i++ ) { cout << " " << " " << setw(14) << w[i] << " " << setw(14) << x[i] << "\n"; } delete [] w; delete [] x; } cout.precision ( prec ); return; } //****************************************************************************80 void gegenbauer_to_monomial_matrix_test ( ) //****************************************************************************80 // // Purpose: // // gegenbauer_to_monomial_matrix_test() tests gegenbauer_to_monomial_matrix(). // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 16 April 2024 // // Author: // // John Burkardt // { double *A; double alpha; double *g; int i; int j; int n = 5; double *mono; cout << "\n"; cout << "gegenbauer_to_monomial_matrix_test():\n"; cout << " gegenbauer_to_monomial_matrix() evaluates the matrix\n"; cout << " which converts Gegenbauer polyjomial coefficients\n"; cout << " to monomial coefficients.\n"; alpha = 0.5; A = gegenbauer_to_monomial_matrix ( n, alpha ); r8mat_print ( n, n, A, " Gegenbauer to Monomial matrix G:" ); g = new double[n]; for ( i = 0; i < n; i++ ) { for ( j = 0; j < n; j++ ) { g[j] = 0.0; } g[i] = 1.0; mono = r8mat_mv_new ( n, n, A, g ); cout << "\n"; cout << " Monomial form of Gegenbauer polynomial " << i << "\n"; r8poly_print ( i, mono, "" ); delete [] mono; } delete [] A; delete [] g; return; } //****************************************************************************80 void imtqlx_test ( ) //****************************************************************************80 // // Purpose: // // imtqlx_test() tests imtqlx(). // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 10 June 2015 // // Author: // // John Burkardt. // { double angle; double d[5]; double e[5]; int i; double lam[5]; double lam2[5]; int n = 5; double qtz[5]; double r8_pi = 3.141592653589793; double z[5]; cout << "\n"; cout << "IMTQLX_TEST\n"; cout << " IMTQLX takes a symmetric tridiagonal matrix A\n"; cout << " and computes its eigenvalues LAM.\n"; cout << " It also accepts a vector Z and computes Q'*Z,\n"; cout << " where Q is the matrix that diagonalizes A.\n"; for ( i = 0; i < n; i++ ) { d[i] = 2.0; } for ( i = 0; i < n - 1; i++ ) { e[i] = -1.0; } e[n-1] = 0.0; for ( i = 0; i < n; i++ ) { z[i] = 1.0; } // // On input, LAM is D, and QTZ is Z. // for ( i = 0; i < n; i++ ) { lam[i] = d[i]; } for ( i = 0; i < n; i++ ) { qtz[i] = z[i]; } imtqlx ( n, lam, e, qtz ); r8vec_print ( n, lam, " Computed eigenvalues:" ); for ( i = 0; i < n; i++ ) { angle = ( double ) ( i + 1 ) * r8_pi / ( double ) ( 2 * ( n + 1 ) ); lam2[i] = 4.0 * pow ( sin ( angle ), 2 ); } r8vec_print ( n, lam2, " Exact eigenvalues:" ); r8vec_print ( n, z, " Vector Z:" ); r8vec_print ( n, qtz, " Vector Q'*Z:" ); return; } //****************************************************************************80 void monomial_to_gegenbauer_matrix_test ( ) //****************************************************************************80 // // Purpose: // // monomial_to_gegenbauer_matrix_test() tests monomial_to_gegenbauer_matrix(). // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 24 April 2024 // // Author: // // John Burkardt // { double alpha; double *G; double *Ginv; double *I; int n = 5; cout << "\n"; cout << "monomial_to_gegenbauer_matrix_test():\n"; cout << " monomial_to_gegenbauer_matrix() evaluates the matrix\n"; cout << " which converts monomial polynomial coefficients\n"; cout << " to Gegenbauer coefficients.\n"; alpha = 0.5; cout << " Using parameter alpha = " << alpha << "\n"; G = gegenbauer_to_monomial_matrix ( n, alpha ); r8mat_print ( n, n, G, " Gegenbauer to Monomial matrix G:" ); Ginv = monomial_to_gegenbauer_matrix ( n, alpha ); r8mat_print ( n, n, Ginv, " Monomial to Gegenbauer matrix Ginv:" ); I = r8mat_mm_new ( n, n, n, G, Ginv ); r8mat_print ( n, n, I, " I = G * Ginv:" ); delete [] G; delete [] Ginv; delete [] I; return; } //****************************************************************************80 void r8_hyper_2f1_test ( ) //****************************************************************************80 // // Purpose: // // r8_hyper_2f1_test() tests r8_hyper_2f1(). // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 09 February 2008 // // Author: // // John Burkardt // { double a; double b; double c; double fx; double fx2; int n_data; double x; cout << "\n"; cout << " R8_HYPER_2F1_TEST:\n"; cout << " R8_HYPER_2F1 evaluates the hypergeometric function 2F1.\n"; cout << "\n"; cout << " A B C X "; cout << " 2F1 2F1 DIFF\n"; cout << " "; cout << "(tabulated) (computed)\n"; cout << "\n"; n_data = 0; for ( ; ; ) { hyper_2f1_values ( n_data, a, b, c, x, fx ); if ( n_data == 0 ) { break; } fx2 = r8_hyper_2f1 ( a, b, c, x ); cout << " " << setw(6) << setprecision(2) << a << " " << setw(6) << setprecision(2) << b << " " << setw(6) << setprecision(2) << c << " " << setw(6) << setprecision(2) << x << " " << setw(24) << setprecision(16) << fx << " " << setw(24) << setprecision(16) << fx2 << " " << setw(10) << setprecision(4) << fabs ( fx - fx2 ) << "\n"; } return; } //****************************************************************************80 void r8_uniform_ab_test ( ) //****************************************************************************80 // // Purpose: // // r8_uniform_ab_test() tests r8_uniform_ab(). // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 11 April 2007 // // Author: // // John Burkardt // { double a; double b; double c; int i; int seed; b = 10.0; c = 25.0; seed = 17; cout << "\n"; cout << "r8_uniform_ab_test():\n"; cout << " r8_uniform_ab() produces a random real in a given range.\n"; cout << "\n"; cout << " Using range " << b << " <= A <= " << c << ".\n"; cout << "\n"; cout << "\n"; cout << " I A\n"; cout << "\n"; for ( i = 0; i < 10; i++ ) { a = r8_uniform_ab ( b, c, seed ); cout << setw ( 6 ) << i << " " << setw ( 10 ) << a << "\n"; } return; } //****************************************************************************80 double *r8mat_mm_new ( int n1, int n2, int n3, double a[], double b[] ) //****************************************************************************80 // // Purpose: // // r8mat_mm_new() multiplies two matrices. // // Discussion: // // An R8MAT is a doubly dimensioned array of R8 values, stored as a vector // in column-major order. // // For this routine, the result is returned as the function value. // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 18 October 2005 // // Author: // // John Burkardt // // Parameters: // // Input, int N1, N2, N3, the order of the matrices. // // Input, double A[N1*N2], double B[N2*N3], the matrices to multiply. // // Output, double R8MAT_MM_NEW[N1*N3], the product matrix C = A * B. // { double *c; int i; int j; int k; c = new double[n1*n3]; for ( i = 0; i < n1; i++ ) { for ( j = 0; j < n3; j++ ) { c[i+j*n1] = 0.0; for ( k = 0; k < n2; k++ ) { c[i+j*n1] = c[i+j*n1] + a[i+k*n1] * b[k+j*n2]; } } } return c; } //****************************************************************************80 double *r8mat_mv_new ( int m, int n, double a[], double x[] ) //****************************************************************************80 // // Purpose: // // r8mat_mv_new() multiplies a matrix times a vector. // // Discussion: // // An R8MAT is a doubly dimensioned array of R8 values, stored as a vector // in column-major order. // // For this routine, the result is returned as the function value. // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 11 April 2007 // // Author: // // John Burkardt // // Input: // // int M, N, the number of rows and columns of the matrix. // // double A[M,N], the M by N matrix. // // double X[N], the vector to be multiplied by A. // // Output: // // double R8MAT_MV_NEW[M], the product A*X. // { int i; int j; double *y; y = new double[m]; for ( i = 0; i < m; i++ ) { y[i] = 0.0; for ( j = 0; j < n; j++ ) { y[i] = y[i] + a[i+j*m] * x[j]; } } return y; } //****************************************************************************80 void r8mat_print ( int m, int n, double a[], string title ) //****************************************************************************80 // // Purpose: // // r8mat_print() prints an R8MAT. // // Discussion: // // An R8MAT is a doubly dimensioned array of R8 values, stored as a vector // in column-major order. // // Entry A(I,J) is stored as A[I+J*M] // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 10 September 2009 // // Author: // // John Burkardt // // Input: // // int M, the number of rows in A. // // int N, the number of columns in A. // // double A[M*N], the M by N matrix. // // string TITLE, a title. // { r8mat_print_some ( m, n, a, 1, 1, m, n, title ); return; } //****************************************************************************80 void r8mat_print_some ( int m, int n, double a[], int ilo, int jlo, int ihi, int jhi, string title ) //****************************************************************************80 // // Purpose: // // r8mat_print_some() prints some of an R8MAT. // // Discussion: // // An R8MAT is a doubly dimensioned array of R8 values, stored as a vector // in column-major order. // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 26 June 2013 // // Author: // // John Burkardt // // Input: // // int M, the number of rows of the matrix. // M must be positive. // // int N, the number of columns of the matrix. // N must be positive. // // double A[M*N], the matrix. // // int ILO, JLO, IHI, JHI, designate the first row and // column, and the last row and column to be printed. // // string TITLE, a title. // { # define INCX 5 int i; int i2hi; int i2lo; int j; int j2hi; int j2lo; cout << "\n"; cout << title << "\n"; if ( m <= 0 || n <= 0 ) { cout << "\n"; cout << " (None)\n"; return; } // // Print the columns of the matrix, in strips of 5. // for ( j2lo = jlo; j2lo <= jhi; j2lo = j2lo + INCX ) { j2hi = j2lo + INCX - 1; if ( n < j2hi ) { j2hi = n; } if ( jhi < j2hi ) { j2hi = jhi; } cout << "\n"; // // For each column J in the current range... // // Write the header. // cout << " Col: "; for ( j = j2lo; j <= j2hi; j++ ) { cout << setw(7) << j - 1 << " "; } cout << "\n"; cout << " Row\n"; cout << "\n"; // // Determine the range of the rows in this strip. // if ( 1 < ilo ) { i2lo = ilo; } else { i2lo = 1; } if ( ihi < m ) { i2hi = ihi; } else { i2hi = m; } for ( i = i2lo; i <= i2hi; i++ ) { // // Print out (up to) 5 entries in row I, that lie in the current strip. // cout << setw(5) << i - 1 << ": "; for ( j = j2lo; j <= j2hi; j++ ) { cout << setw(12) << a[i-1+(j-1)*m] << " "; } cout << "\n"; } } return; # undef INCX } //****************************************************************************80 void r8poly_print ( int n, double a[], string title ) //****************************************************************************80 // // Purpose: // // r8poly_print() prints a polynomial. // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 15 July 2015 // // Author: // // John Burkardt // // Input: // // int N, the dimension of A. // // double A[N+1], the polynomial coefficients. // A(0) is the constant term and // A(N) is the coefficient of X^N. // // string TITLE, a title. // { int i; double mag; char plus_minus; if ( 0 < title.length ( ) ) { cout << "\n"; cout << title << "\n"; } cout << "\n"; if ( n < 0 ) { cout << " p(x) = 0\n"; return; } if ( a[n] < 0.0 ) { plus_minus = '-'; } else { plus_minus = ' '; } mag = fabs ( a[n] ); if ( 2 <= n ) { cout << " p(x) = " << plus_minus << setw(14) << mag << " * x ^ " << n << "\n"; } else if ( n == 1 ) { cout << " p(x) = " << plus_minus << setw(14) << mag << " * x\n"; } else if ( n == 0 ) { cout << " p(x) = " << plus_minus << setw(14) << mag << "\n"; } for ( i = n - 1; 0 <= i; i-- ) { if ( a[i] < 0.0 ) { plus_minus = '-'; } else { plus_minus = '+'; } mag = fabs ( a[i] ); if ( mag != 0.0 ) { if ( 2 <= i ) { cout << " " << plus_minus << setw(14) << mag << " * x ^ " << i << "\n"; } else if ( i == 1 ) { cout << " " << plus_minus << setw(14) << mag << " * x\n"; } else if ( i == 0 ) { cout << " " << plus_minus << setw(14) << mag << "\n"; } } } return; } //****************************************************************************80 void timestamp ( ) //****************************************************************************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 MIT license. // // Modified: // // 19 March 2018 // // Author: // // John Burkardt // { # define TIME_SIZE 40 static char time_buffer[TIME_SIZE]; const struct std::tm *tm_ptr; std::time_t now; now = std::time ( NULL ); tm_ptr = std::localtime ( &now ); std::strftime ( time_buffer, TIME_SIZE, "%d %B %Y %I:%M:%S %p", tm_ptr ); std::cout << time_buffer << "\n"; return; # undef TIME_SIZE }