# include # include # include # include # include using namespace std; # include "hypersphere_monte_carlo.hpp" //****************************************************************************80 double hypersphere01_area ( int m ) //****************************************************************************80 // // Purpose: // // HYPERSPHERE01_AREA returns the area of the unit hypersphere. // // Discussion: // // M Area // // 2 2 * PI // 3 4 * PI // 4 ( 2 / 1) * PI^2 // 5 ( 8 / 3) * PI^2 // 6 ( 1 / 1) * PI^3 // 7 (16 / 15) * PI^3 // 8 ( 1 / 3) * PI^4 // 9 (32 / 105) * PI^4 // 10 ( 1 / 12) * PI^5 // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 04 January 2014 // // Author: // // John Burkardt // // Parameters: // // Input, int M, the spatial dimension. // // Output, double HYPERSPHERE01_AREA, the area. // { double area; int i; int m_half; const double r8_pi = 3.141592653589793; if ( ( m % 2 ) == 0 ) { m_half = m / 2; area = 2.0 * pow ( r8_pi, m_half ); for ( i = 1; i <= m_half - 1; i++ ) { area = area / ( double ) ( i ); } } else { m_half = ( m - 1 ) / 2; area = pow ( r8_pi, m_half ) * pow ( 2.0, m ); for ( i = m_half + 1; i <= 2 * m_half; i++ ) { area = area / ( double ) ( i ); } } return area; } //****************************************************************************80 double hypersphere01_monomial_integral ( int m, int e[] ) //****************************************************************************80 // // Purpose: // // HYPERSPHERE01_MONOMIAL_INTEGRAL: monomial integrals on unit hypersphere. // // Discussion: // // The integration region is // // sum ( 1 <= I <= M ) X(I)^2 = 1. // // The monomial is F(X) = product ( 1 <= I <= M ) X(I)^E(I). // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 04 January 2014 // // Author: // // John Burkardt // // Reference: // // Philip Davis, Philip Rabinowitz, // Methods of Numerical Integration, // Second Edition, // Academic Press, 1984, page 263. // // Parameters: // // Input, int M, the spatial dimension. // // Input, int E[M], the exponents. Each exponent must be nonnegative. // // Output, double HYPERSPHERE01_MONOMIAL_INTEGRAL, the integral. // { double arg; int i; double integral; for ( i = 0; i < m; i++ ) { if ( e[i] < 0 ) { cout << "\n"; cout << "HYPERSPHERE01_MONOMIAL_INTEGRAL - Fatal error!\n"; cout << " All exponents must be nonnegative.\n"; cout << " E[" << i << "] = " << e[0] << "\n"; exit ( 1 ); } } for ( i = 0; i < m; i++ ) { if ( ( e[i] % 2 ) == 1 ) { integral = 0.0; return integral; } } integral = 2.0; for ( i = 0; i < m; i++ ) { arg = 0.5 * ( double ) ( e[i] + 1 ); integral = integral * tgamma ( arg ); } arg = 0.5 * ( double ) ( i4vec_sum ( m, e ) + m ); integral = integral / tgamma ( arg ); return integral; } //****************************************************************************80 double *hypersphere01_sample ( int m, int n, int &seed ) //****************************************************************************80 // // Purpose: // // HYPERSPHERE01_SAMPLE uniformly samples the surface of the unit hypersphere. // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 04 January 2014 // // Author: // // John Burkardt // // Reference: // // Russell Cheng, // Random Variate Generation, // in Handbook of Simulation, // edited by Jerry Banks, // Wiley, 1998, pages 168. // // Reuven Rubinstein, // Monte Carlo Optimization, Simulation, and Sensitivity // of Queueing Networks, // Krieger, 1992, // ISBN: 0894647644, // LC: QA298.R79. // // Parameters: // // Input, int M, the spatial dimension. // // Input, int N, the number of points. // // Input/output, int &SEED, a seed for the random // number generator. // // Output, double X[M*N], the points. // { int i; int j; double norm; double *x; x = r8mat_normal_01_new ( m, n, seed ); for ( j = 0; j < n; j++ ) { // // Compute the length of the vector. // norm = 0.0; for ( i = 0; i < m; i++ ) { norm = norm + pow ( x[i+j*m], 2 ); } norm = sqrt ( norm ); // // Normalize the vector. // for ( i = 0; i < m; i++ ) { x[i+j*m] = x[i+j*m] / norm; } } return x; } //****************************************************************************80 int i4vec_sum ( int n, int a[] ) //****************************************************************************80 // // Purpose: // // I4VEC_SUM sums the entries of an I4VEC. // // Discussion: // // An I4VEC is a vector of I4's. // // Example: // // Input: // // A = ( 1, 2, 3, 4 ) // // Output: // // I4VEC_SUM = 10 // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 26 May 1999 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of entries in the vector. // // Input, int A[N], the vector to be summed. // // Output, int I4VEC_SUM, the sum of the entries of A. // { int i; int sum; sum = 0; for ( i = 0; i < n; i++ ) { sum = sum + a[i]; } return sum; } //****************************************************************************80 double *monomial_value ( int m, int n, int e[], double x[] ) //****************************************************************************80 // // Purpose: // // MONOMIAL_VALUE evaluates a monomial. // // Discussion: // // This routine evaluates a monomial of the form // // product ( 1 <= i <= m ) x(i)^e(i) // // where the exponents are nonnegative integers. Note that // if the combination 0^0 is encountered, it should be treated // as 1. // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 05 May 2007 // // Author: // // John Burkardt // // Parameters: // // Input, int M, the spatial dimension. // // Input, int POINT_NUM, the number of points at which the // monomial is to be evaluated. // // Input, int E[M], the exponents. // // Input, double X[M*N], the point coordinates. // // Output, double MONOMIAL_VALUE[N], the value of the monomial. // { int i; int j; double *v; v = new double[n]; for ( j = 0; j < n; j++ ) { v[j] = 1.0; } for ( i = 0; i < m; i++ ) { if ( 0 != e[i] ) { for ( j = 0; j < n; j++ ) { v[j] = v[j] * pow ( x[i+j*m], e[i] ); } } } return v; } //****************************************************************************80 double r8_uniform_01 ( int &seed ) //****************************************************************************80 // // Purpose: // // R8_UNIFORM_01 returns a unit pseudorandom R8. // // Discussion: // // This routine implements the recursion // // seed = ( 16807 * seed ) mod ( 2^31 - 1 ) // u = seed / ( 2^31 - 1 ) // // The integer arithmetic never requires more than 32 bits, // including a sign bit. // // If the initial seed is 12345, then the first three computations are // // Input Output R8_UNIFORM_01 // SEED SEED // // 12345 207482415 0.096616 // 207482415 1790989824 0.833995 // 1790989824 2035175616 0.947702 // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 09 April 2012 // // Author: // // John Burkardt // // Reference: // // Paul Bratley, Bennett Fox, Linus Schrage, // A Guide to Simulation, // Second Edition, // Springer, 1987, // ISBN: 0387964673, // LC: QA76.9.C65.B73. // // Bennett Fox, // Algorithm 647: // Implementation and Relative Efficiency of Quasirandom // Sequence Generators, // ACM Transactions on Mathematical Software, // Volume 12, Number 4, December 1986, pages 362-376. // // Pierre L'Ecuyer, // Random Number Generation, // in Handbook of Simulation, // edited by Jerry Banks, // Wiley, 1998, // ISBN: 0471134031, // LC: T57.62.H37. // // Peter Lewis, Allen Goodman, James Miller, // A Pseudo-Random Number Generator for the System/360, // IBM Systems Journal, // Volume 8, Number 2, 1969, pages 136-143. // // Parameters: // // Input/output, int &SEED, the "seed" value. Normally, this // value should not be 0. On output, SEED has been updated. // // Output, double R8_UNIFORM_01, a new pseudorandom variate, // strictly between 0 and 1. // { int i4_huge = 2147483647; int k; double r; if ( seed == 0 ) { cerr << "\n"; cerr << "R8_UNIFORM_01 - Fatal error!\n"; cerr << " Input value of SEED = 0.\n"; exit ( 1 ); } k = seed / 127773; seed = 16807 * ( seed - k * 127773 ) - k * 2836; if ( seed < 0 ) { seed = seed + i4_huge; } r = ( double ) ( seed ) * 4.656612875E-10; return r; } //****************************************************************************80 double *r8mat_normal_01_new ( int m, int n, int &seed ) //****************************************************************************80 // // Purpose: // // R8MAT_NORMAL_01_NEW returns a unit pseudonormal R8MAT. // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 09 December 2009 // // Author: // // John Burkardt // // Reference: // // Paul Bratley, Bennett Fox, Linus Schrage, // A Guide to Simulation, // Springer Verlag, pages 201-202, 1983. // // Bennett Fox, // Algorithm 647: // Implementation and Relative Efficiency of Quasirandom // Sequence Generators, // ACM Transactions on Mathematical Software, // Volume 12, Number 4, pages 362-376, 1986. // // Peter Lewis, Allen Goodman, James Miller, // A Pseudo-Random Number Generator for the System/360, // IBM Systems Journal, // Volume 8, pages 136-143, 1969. // // Parameters: // // Input, int M, N, the number of rows and columns in the array. // // Input/output, int &SEED, the "seed" value, which should NOT be 0. // On output, SEED has been updated. // // Output, double R8MAT_NORMAL_01_NEW[M*N], the array of pseudonormal values. // { double *r; r = r8vec_normal_01_new ( m * n, seed ); return r; } //****************************************************************************80 double *r8vec_normal_01_new ( int n, int &seed ) //****************************************************************************80 // // Purpose: // // R8VEC_NORMAL_01_NEW returns a unit pseudonormal R8VEC. // // Discussion: // // An R8VEC is a vector of R8's. // // The standard normal probability distribution function (PDF) has // mean 0 and standard deviation 1. // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 06 August 2013 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of values desired. // // Input/output, int &SEED, a seed for the random number generator. // // Output, double R8VEC_NORMAL_01_NEW[N], a sample of the standard normal PDF. // // Local parameters: // // Local, double R[N+1], is used to store some uniform random values. // Its dimension is N+1, but really it is only needed to be the // smallest even number greater than or equal to N. // // Local, int X_LO, X_HI, records the range of entries of // X that we need to compute. // { int i; int m; const double pi = 3.141592653589793; double *r; double *x; int x_hi; int x_lo; x = new double[n]; // // Record the range of X we need to fill in. // x_lo = 1; x_hi = n; // // If we need just one new value, do that here to avoid null arrays. // if ( x_hi - x_lo + 1 == 1 ) { r = r8vec_uniform_01_new ( 2, seed ); x[x_hi-1] = sqrt ( -2.0 * log ( r[0] ) ) * cos ( 2.0 * pi * r[1] ); delete [] r; } // // If we require an even number of values, that's easy. // else if ( ( x_hi - x_lo + 1 ) % 2 == 0 ) { m = ( x_hi - x_lo + 1 ) / 2; r = r8vec_uniform_01_new ( 2*m, seed ); for ( i = 0; i <= 2*m-2; i = i + 2 ) { x[x_lo+i-1] = sqrt ( -2.0 * log ( r[i] ) ) * cos ( 2.0 * pi * r[i+1] ); x[x_lo+i ] = sqrt ( -2.0 * log ( r[i] ) ) * sin ( 2.0 * pi * r[i+1] ); } delete [] r; } // // If we require an odd number of values, we generate an even number, // and handle the last pair specially, storing one in X(N), and // saving the other for later. // else { x_hi = x_hi - 1; m = ( x_hi - x_lo + 1 ) / 2 + 1; r = r8vec_uniform_01_new ( 2*m, seed ); for ( i = 0; i <= 2*m-4; i = i + 2 ) { x[x_lo+i-1] = sqrt ( -2.0 * log ( r[i] ) ) * cos ( 2.0 * pi * r[i+1] ); x[x_lo+i ] = sqrt ( -2.0 * log ( r[i] ) ) * sin ( 2.0 * pi * r[i+1] ); } i = 2*m - 2; x[x_lo+i-1] = sqrt ( -2.0 * log ( r[i] ) ) * cos ( 2.0 * pi * r[i+1] ); delete [] r; } return x; } //****************************************************************************80 double r8vec_sum ( int n, double a[] ) //****************************************************************************80 // // Purpose: // // R8VEC_SUM returns the sum of an R8VEC. // // Discussion: // // An R8VEC is a vector of R8's. // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 15 October 2004 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of entries in the vector. // // Input, double A[N], the vector. // // Output, double R8VEC_SUM, the sum of the vector. // { int i; double value; value = 0.0; for ( i = 0; i < n; i++ ) { value = value + a[i]; } return value; } //****************************************************************************80 double *r8vec_uniform_01_new ( int n, int &seed ) //****************************************************************************80 // // Purpose: // // R8VEC_UNIFORM_01_NEW returns a new unit pseudorandom R8VEC. // // Discussion: // // This routine implements the recursion // // seed = ( 16807 * seed ) mod ( 2^31 - 1 ) // u = seed / ( 2^31 - 1 ) // // The integer arithmetic never requires more than 32 bits, // including a sign bit. // // Licensing: // // This code is distributed under the MIT license. // // Modified: // // 19 August 2004 // // Author: // // John Burkardt // // Reference: // // Paul Bratley, Bennett Fox, Linus Schrage, // A Guide to Simulation, // Second Edition, // Springer, 1987, // ISBN: 0387964673, // LC: QA76.9.C65.B73. // // Bennett Fox, // Algorithm 647: // Implementation and Relative Efficiency of Quasirandom // Sequence Generators, // ACM Transactions on Mathematical Software, // Volume 12, Number 4, December 1986, pages 362-376. // // Pierre L'Ecuyer, // Random Number Generation, // in Handbook of Simulation, // edited by Jerry Banks, // Wiley, 1998, // ISBN: 0471134031, // LC: T57.62.H37. // // Peter Lewis, Allen Goodman, James Miller, // A Pseudo-Random Number Generator for the System/360, // IBM Systems Journal, // Volume 8, Number 2, 1969, pages 136-143. // // Parameters: // // Input, int N, the number of entries in the vector. // // Input/output, int &SEED, a seed for the random number generator. // // Output, double R8VEC_UNIFORM_01_NEW[N], the vector of pseudorandom values. // { int i; int i4_huge = 2147483647; int k; double *r; if ( seed == 0 ) { cerr << "\n"; cerr << "R8VEC_UNIFORM_01_NEW - Fatal error!\n"; cerr << " Input value of SEED = 0.\n"; exit ( 1 ); } r = new double[n]; for ( i = 0; i < n; i++ ) { k = seed / 127773; seed = 16807 * ( seed - k * 127773 ) - k * 2836; if ( seed < 0 ) { seed = seed + i4_huge; } r[i] = ( double ) ( seed ) * 4.656612875E-10; } return r; } //****************************************************************************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: // // 08 July 2009 // // Author: // // John Burkardt // // Parameters: // // None // { # 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 }