# include # include # include # include # include using namespace std; # include "haar.hpp" //****************************************************************************80 void haar_1d ( int n, double x[] ) //****************************************************************************80 // // Purpose: // // HAAR_1D computes the Haar transform of a vector. // // Discussion: // // For the classical Haar transform, N should be a power of 2. // However, this is not required here. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 06 March 2014 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the dimension of the vector. // // Input/output, double X[N], on input, the vector to be transformed. // On output, the transformed vector. // { int i; int k; double s; double *y; s = sqrt ( 2.0 ); y = new double[n]; // // Initialize. // for ( i = 0; i < n; i++ ) { y[i] = 0.0; } // // Determine K, the largest power of 2 such that K <= N. // k = 1; while ( k * 2 <= n ) { k = k * 2; } while ( 1 < k ) { k = k / 2; for ( i = 0; i < k; i++ ) { y[i] = ( x[2*i] + x[2*i+1] ) / s; y[i+k] = ( x[2*i] - x[2*i+1] ) / s; } for ( i = 0; i < k * 2; i++ ) { x[i] = y[i]; } } delete [] y; return; } //****************************************************************************80 void haar_1d_inverse ( int n, double x[] ) //****************************************************************************80 // // Purpose: // // HAAR_1D_INVERSE computes the inverse Haar transform of a vector. // // Discussion: // // For the classical Haar transform, N should be a power of 2. // However, this is not required here. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 06 March 2014 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the dimension of the vector. // // Input/output, double X[N], on input, the vector to be transformed. // On output, the transformed vector. // { int i; int k; double s; double *y; s = sqrt ( 2.0 ); y = new double[n]; // // Initialize. // for ( i = 0; i < n; i++ ) { y[i] = 0.0; } k = 1; while ( k * 2 <= n ) { for ( i = 0; i < k; i++ ) { y[2*i] = ( x[i] + x[i+k] ) / s; y[2*i+1] = ( x[i] - x[i+k] ) / s; } for ( i = 0; i < k * 2; i++ ) { x[i] = y[i]; } k = k * 2; } delete [] y; return; } //****************************************************************************80 void haar_2d ( int m, int n, double u[] ) //****************************************************************************80 // // Purpose: // // HAAR_2D computes the Haar transform of an array. // // Discussion: // // For the classical Haar transform, M and N should be a power of 2. // However, this is not required here. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 06 March 2014 // // Author: // // John Burkardt // // Parameters: // // Input, int M, N, the dimensions of the array. // // Input/output, double U[M*N], the array to be transformed. // { int i; int j; int k; double s; double *v; s = sqrt ( 2.0 ); v = new double[m*n]; for ( j = 0; j < n; j++ ) { for ( i = 0; i < m; i++ ) { v[i+j*m] = u[i+j*m]; } } // // Determine K, the largest power of 2 such that K <= M. // k = 1; while ( k * 2 <= m ) { k = k * 2; } // // Transform all columns. // while ( 1 < k ) { k = k / 2; for ( j = 0; j < n; j++ ) { for ( i = 0; i < k; i++ ) { v[i +j*m] = ( u[2*i+j*m] + u[2*i+1+j*m] ) / s; v[k+i+j*m] = ( u[2*i+j*m] - u[2*i+1+j*m] ) / s; } } for ( j = 0; j < n; j++ ) { for ( i = 0; i < 2 * k; i++ ) { u[i+j*m] = v[i+j*m]; } } } // // Determine K, the largest power of 2 such that K <= N. // k = 1; while ( k * 2 <= n ) { k = k * 2; } // // Transform all rows. // while ( 1 < k ) { k = k / 2; for ( j = 0; j < k; j++ ) { for ( i = 0; i < m; i++ ) { v[i+( j)*m] = ( u[i+2*j*m] + u[i+(2*j+1)*m] ) / s; v[i+(k+j)*m] = ( u[i+2*j*m] - u[i+(2*j+1)*m] ) / s; } } for ( j = 0; j < 2 * k; j++ ) { for ( i = 0; i < m; i++ ) { u[i+j*m] = v[i+j*m]; } } } delete [] v; return; } //****************************************************************************80 void haar_2d_inverse ( int m, int n, double u[] ) //****************************************************************************80 // // Purpose: // // HAAR_2D_INVERSE inverts the Haar transform of an array. // // Discussion: // // For the classical Haar transform, M and N should be a power of 2. // However, this is not required here. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 06 March 2014 // // Author: // // John Burkardt // // Parameters: // // Input, int M, N, the dimensions of the array. // // Input/output, double U[M*N], the array to be transformed. // { int i; int j; int k; double s; double *v; s = sqrt ( 2.0 ); v = new double[m*n]; for ( j = 0; j < n; j++ ) { for ( i = 0; i < m; i++ ) { v[i+j*m] = u[i+j*m]; } } // // Inverse transform of all rows. // k = 1; while ( k * 2 <= n ) { for ( j = 0; j < k; j++ ) { for ( i = 0; i < m; i++ ) { v[i+(2*j )*m] = ( u[i+j*m] + u[i+(k+j)*m] ) / s; v[i+(2*j+1)*m] = ( u[i+j*m] - u[i+(k+j)*m] ) / s; } } for ( j = 0; j < 2 * k; j++ ) { for ( i = 0; i < m; i++ ) { u[i+j*m] = v[i+j*m]; } } k = k * 2; } // // Inverse transform of all columns. // k = 1; while ( k * 2 <= m ) { for ( j = 0; j < n; j++ ) { for ( i = 0; i < k; i++ ) { v[2*i +j*m] = ( u[i+j*m] + u[k+i+j*m] ) / s; v[2*i+1+j*m] = ( u[i+j*m] - u[k+i+j*m] ) / s; } } for ( j = 0; j < n; j++ ) { for ( i = 0; i < 2 * k; i++ ) { u[i+j*m] = v[i+j*m]; } } k = k * 2; } delete [] v; return; } //****************************************************************************80 int i4_max ( int i1, int i2 ) //****************************************************************************80 // // Purpose: // // I4_MAX returns the maximum of two I4's. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 13 October 1998 // // Author: // // John Burkardt // // Parameters: // // Input, int I1, I2, are two integers to be compared. // // Output, int I4_MAX, the larger of I1 and I2. // { int value; if ( i2 < i1 ) { value = i1; } else { value = i2; } return value; } //****************************************************************************80 int i4_min ( int i1, int i2 ) //****************************************************************************80 // // Purpose: // // I4_MIN returns the minimum of two I4's. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 13 October 1998 // // Author: // // John Burkardt // // Parameters: // // Input, int I1, I2, two integers to be compared. // // Output, int I4_MIN, the smaller of I1 and I2. // { int value; if ( i1 < i2 ) { value = i1; } else { value = i2; } return value; } //****************************************************************************80 double *r8mat_copy_new ( int m, int n, double a1[] ) //****************************************************************************80 // // Purpose: // // R8MAT_COPY_NEW copies one R8MAT to a "new" R8MAT. // // Discussion: // // An R8MAT is a doubly dimensioned array of R8's, which // may be stored as a vector in column-major order. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 03 July 2008 // // Author: // // John Burkardt // // Parameters: // // Input, int M, N, the number of rows and columns. // // Input, double A1[M*N], the matrix to be copied. // // Output, double R8MAT_COPY_NEW[M*N], the copy of A1. // { double *a2; int i; int j; a2 = new double[m*n]; for ( j = 0; j < n; j++ ) { for ( i = 0; i < m; i++ ) { a2[i+j*m] = a1[i+j*m]; } } return a2; } //****************************************************************************80 double r8mat_dif_fro ( int m, int n, double a[], double b[] ) //****************************************************************************80 // // Purpose: // // R8MAT_DIF_FRO returns the Frobenius norm of the difference of R8MAT's. // // Discussion: // // An R8MAT is a doubly dimensioned array of double precision values, which // may be stored as a vector in column-major order. // // The Frobenius norm is defined as // // R8MAT_NORM_FRO = sqrt ( // sum ( 1 <= I <= M ) sum ( 1 <= j <= N ) A(I,J)^2 ) // // The matrix Frobenius norm is not derived from a vector norm, but // is compatible with the vector L2 norm, so that: // // r8vec_norm_l2 ( A * x ) <= r8mat_norm_fro ( A ) * r8vec_norm_l2 ( x ). // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 14 September 2006 // // Author: // // John Burkardt // // Parameters: // // Input, int M, the number of rows in A. // // Input, int N, the number of columns in A. // // Input, double A[M*N], double B[M*N], the matrices for which we // want the Frobenius norm of the difference. // // Output, double R8MAT_DIF_FRO, the Frobenius norm of ( A - B ). // { int i; int j; double value; value = 0.0; for ( j = 0; j < n; j++ ) { for ( i = 0; i < m; i++ ) { value = value + pow ( a[i+j*m] - b[i+j*m], 2 ); } } value = sqrt ( value ); return value; } //****************************************************************************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 GNU LGPL license. // // Modified: // // 10 September 2009 // // Author: // // John Burkardt // // Parameters: // // Input, int M, the number of rows in A. // // Input, int N, the number of columns in A. // // Input, double A[M*N], the M by N matrix. // // Input, 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 GNU LGPL license. // // Modified: // // 20 August 2010 // // Author: // // John Burkardt // // Parameters: // // Input, int M, the number of rows of the matrix. // M must be positive. // // Input, int N, the number of columns of the matrix. // N must be positive. // // Input, double A[M*N], the matrix. // // Input, int ILO, JLO, IHI, JHI, designate the first row and // column, and the last row and column to be printed. // // Input, 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; j2hi = i4_min ( j2hi, n ); j2hi = i4_min ( 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. // i2lo = i4_max ( ilo, 1 ); i2hi = i4_min ( ihi, 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 double *r8mat_uniform_01_new ( int m, int n, int &seed ) //****************************************************************************80 // // Purpose: // // R8MAT_UNIFORM_01_NEW returns a unit pseudorandom R8MAT. // // Discussion: // // An R8MAT is a doubly dimensioned array of R8's, stored as a vector // in column-major order. // // This routine implements the recursion // // seed = 16807 * seed mod ( 2^31 - 1 ) // unif = seed / ( 2^31 - 1 ) // // The integer arithmetic never requires more than 32 bits, // including a sign bit. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 03 October 2005 // // 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. // // Philip 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. // // Input/output, int &SEED, the "seed" value. Normally, this // value should not be 0, otherwise the output value of SEED // will still be 0, and R8_UNIFORM will be 0. On output, SEED has // been updated. // // Output, double R8MAT_UNIFORM_01_NEW[M*N], a matrix of pseudorandom values. // { int i; int j; int k; double *r; r = new double[m*n]; for ( j = 0; j < n; j++ ) { for ( i = 0; i < m; i++ ) { k = seed / 127773; seed = 16807 * ( seed - k * 127773 ) - k * 2836; if ( seed < 0 ) { seed = seed + 2147483647; } r[i+j*m] = ( double ) ( seed ) * 4.656612875E-10; } } return r; } //****************************************************************************80 double *r8vec_copy_new ( int n, double a1[] ) //****************************************************************************80 // // Purpose: // // R8VEC_COPY_NEW copies an R8VEC to a "new" R8VEC. // // Discussion: // // An R8VEC is a vector of R8's. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 03 July 2008 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of entries in the vectors. // // Input, double A1[N], the vector to be copied. // // Output, double R8VEC_COPY_NEW[N], the copy of A1. // { double *a2; int i; a2 = new double[n]; for ( i = 0; i < n; i++ ) { a2[i] = a1[i]; } return a2; } //****************************************************************************80 double r8vec_diff_norm ( int n, double a[], double b[] ) //****************************************************************************80 // // Purpose: // // R8VEC_DIFF_NORM returns the L2 norm of the difference of R8VEC's. // // Discussion: // // An R8VEC is a vector of R8's. // // The vector L2 norm is defined as: // // R8VEC_NORM_L2 = sqrt ( sum ( 1 <= I <= N ) A(I)^2 ). // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 24 June 2011 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of entries in A. // // Input, double A[N], B[N], the vectors. // // Output, double R8VEC_DIFF_NORM, the L2 norm of A - B. // { int i; double value; value = 0.0; for ( i = 0; i < n; i++ ) { value = value + ( a[i] - b[i] ) * ( a[i] - b[i] ); } value = sqrt ( value ); return value; } //****************************************************************************80 double *r8vec_linspace_new ( int n, double a_first, double a_last ) //****************************************************************************80 // // Purpose: // // R8VEC_LINSPACE_NEW creates a vector of linearly spaced values. // // Discussion: // // An R8VEC is a vector of R8's. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 14 March 2011 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of entries in the vector. // // Input, double A_FIRST, A_LAST, the first and last entries. // // Output, double R8VEC_ONES_NEW[N], a vector of linearly spaced data. // { double *a; int i; a = new double[n]; if ( n == 1 ) { a[0] = ( a_first + a_last ) / 2.0; } else { for ( i = 0; i < n; i++ ) { a[i] = ( ( double ) ( n - 1 - i ) * a_first + ( double ) ( i ) * a_last ) / ( double ) ( n - 1 ); } } return a; } //****************************************************************************80 double *r8vec_ones_new ( int n ) //****************************************************************************80 // // Purpose: // // R8VEC_ONES_NEW creates a vector of 1's. // // Discussion: // // An R8VEC is a vector of R8's. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 14 March 2011 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of entries in the vector. // // Output, double R8VEC_ONES_NEW[N], a vector of 1's. // { double *a; int i; a = new double[n]; for ( i = 0; i < n; i++ ) { a[i] = 1.0; } return a; } //****************************************************************************80 void r8vec_transpose_print ( int n, double a[], string title ) //****************************************************************************80 // // Purpose: // // R8VEC_TRANSPOSE_PRINT prints an R8VEC "transposed". // // Discussion: // // An R8VEC is a vector of R8's. // // Example: // // A = (/ 1.0, 2.1, 3.2, 4.3, 5.4, 6.5, 7.6, 8.7, 9.8, 10.9, 11.0 /) // TITLE = 'My vector: ' // // My vector: // // 1.0 2.1 3.2 4.3 5.4 // 6.5 7.6 8.7 9.8 10.9 // 11.0 // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 12 November 2010 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of components of the vector. // // Input, double A[N], the vector to be printed. // // Input, string TITLE, a title. // { int i; int ihi; int ilo; cout << "\n"; cout << title << "\n"; cout << "\n"; if ( n <= 0 ) { cout << " (Empty)\n"; return; } for ( ilo = 0; ilo < n; ilo = ilo + 5 ) { ihi = i4_min ( ilo + 5, n ); for ( i = ilo; i < ihi; i++ ) { cout << " " << setw(12) << a[i]; } cout << "\n"; } return; } //****************************************************************************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 GNU LGPL 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 GNU LGPL 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; size_t len; std::time_t now; now = std::time ( NULL ); tm_ptr = std::localtime ( &now ); len = std::strftime ( time_buffer, TIME_SIZE, "%d %B %Y %I:%M:%S %p", tm_ptr ); std::cout << time_buffer << "\n"; return; # undef TIME_SIZE }