# include # include # include # include # include "r83v.h" int main ( ); void r83v_cg_test ( ); void r83v_copy_test ( ); void r83v_cr_fa_test ( ); void r83v_cr_sl_test ( ); void r83v_cr_sls_test ( ); void r83v_dif2_test ( ); void r83v_fs_test ( );; void r83v_gs_sl_test ( ); void r83v_indicator_test ( ); void r83v_jac_sl_test ( ); void r83v_mtv_test ( ); void r83v_mv_test ( ); void r83v_print_test ( ); void r83v_print_some_test ( ); void r83v_random_test ( ); void r83v_res_test ( ); void r83v_to_r8ge_test ( ); void r83v_to_r8vec_test ( ); void r83v_transpose_test ( ); void r83v_zeros_test ( ); void r8ge_to_r83v_test ( ); void r8vec_to_r83v_test ( ); void timestamp ( ); /******************************************************************************/ int main ( ) /******************************************************************************/ /* Purpose: r83v_test() tests r83v(). Licensing: This code is distributed under the MIT license. Modified: 16 February 2016 Author: John Burkardt */ { timestamp ( ); printf ( "\n" ); printf ( "R83V_TEST\n" ); printf ( " C version\n" ); printf ( " Test the R83V library.\n" ); r83v_cg_test ( ); r83v_copy_test ( ); r83v_cr_fa_test ( ); r83v_cr_sl_test ( ); r83v_cr_sls_test ( ); r83v_dif2_test ( ); r83v_fs_test ( ); r83v_gs_sl_test ( ); r83v_indicator_test ( ); r83v_jac_sl_test ( ); r83v_mtv_test ( ); r83v_mv_test ( ); r83v_print_test ( ); r83v_print_some_test ( ); r83v_random_test ( ); r83v_res_test ( ); r83v_to_r8ge_test ( ); r83v_to_r8vec_test ( ); r83v_transpose_test ( ); r83v_zeros_test ( ); r8ge_to_r83v_test ( ); r8vec_to_r83v_test ( ); /* Terminate. */ printf ( "\n" ); printf ( "r83v_test():\n" ); printf ( " Normal end of execution.\n" ); printf ( "\n" ); timestamp ( ); return 0; } /******************************************************************************/ void r83v_cg_test ( ) /******************************************************************************/ /* Purpose: R83V_CG_TEST tests R83V_CG. Licensing: This code is distributed under the MIT license. Modified: 14 February 2016 Author: John Burkardt */ { double *a1; double *a2; double *a3; double *b; double e_norm; int i; int n; double *r; double r_norm; int seed; double *x1; double *x2; printf ( "\n" ); printf ( "R83_CG_TEST\n" ); printf ( " R83_CG applies CG to an R83 matrix.\n" ); n = 10; /* Let A be the -1 2 -1 matrix. */ a1 = ( double * ) malloc ( ( n - 1 ) * sizeof ( double ) ); a2 = ( double * ) malloc ( n * sizeof ( double ) ); a3 = ( double * ) malloc ( ( n - 1 ) * sizeof ( double ) ); r83v_dif2 ( n, n, a1, a2, a3 ); /* Choose a random solution. */ seed = 123456789; x1 = r8vec_uniform_01_new ( n, &seed ); /* Compute the corresponding right hand side. */ b = r83v_mv ( n, n, a1, a2, a3, x1 ); /* Call the CG routine. */ x2 = ( double * ) malloc ( n * sizeof ( double ) ); for ( i = 0; i < n; i++ ) { x2[i] = 1.0; } r83v_cg ( n, a1, a2, a3, b, x2 ); /* Compute the residual. */ r = r83v_res ( n, n, a1, a2, a3, x2, b ); r_norm = r8vec_norm ( n, r ); /* Compute the error. */ e_norm = r8vec_norm_affine ( n, x1, x2 ); /* Report. */ printf ( "\n" ); printf ( " Number of variables N = %d\n", n ); printf ( " Norm of residual ||Ax-b|| = %g\n", r_norm ); printf ( " Norm of error ||x1-x2|| = %g\n", e_norm ); /* Free memory. */ free ( a1 ); free ( a2 ); free ( a3 ); free ( b ); free ( r ); free ( x1 ); free ( x2 ); return; } /******************************************************************************/ void r83v_copy_test ( ) /******************************************************************************/ /* Purpose: R83V_COPY_TEST tests R83V_COPY. Licensing: This code is distributed under the MIT license. Modified: 13 February 2016 Author: John Burkardt */ { double *a1; double *a2; double *a3; int ahi; double *b1; double *b2; double *b3; int bhi; int chi; int i; int m; int n; printf ( "\n" ); printf ( "R83V_COPY_TEST\n" ); printf ( " R83V_COPY copies an R83V matrix.\n" ); printf ( " We check three cases, MN.\n" ); for ( i = 1; i <= 3; i++ ) { if ( i == 1 ) { m = 3; n = 5; } else if ( i == 2 ) { m = 5; n = 5; } else if ( i == 3 ) { m = 5; n = 3; } ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a1 = ( double * ) malloc ( ahi * sizeof ( double ) ); a2 = ( double * ) malloc ( bhi * sizeof ( double ) ); a3 = ( double * ) malloc ( chi * sizeof ( double ) ); r83v_indicator ( m, n, a1, a2, a3 ); r83v_print ( m, n, a1, a2, a3, " R83V matrix A:" ); b1 = ( double * ) malloc ( ahi * sizeof ( double ) ); b2 = ( double * ) malloc ( bhi * sizeof ( double ) ); b3 = ( double * ) malloc ( chi * sizeof ( double ) ); r83v_copy ( m, n, a1, a2, a3, b1, b2, b3 ); r83v_print ( m, n, b1, b2, b3, " B = copy of A:" ); free ( a1 ); free ( a2 ); free ( a3 ); free ( b1 ); free ( b2 ); free ( b3 ); } return; } /******************************************************************************/ void r83v_cr_fa_test ( ) /******************************************************************************/ /* Purpose: R83V_CR_FA_TEST tests R83V_CR_FA. Licensing: This code is distributed under the MIT license. Modified: 15 February 2016 Author: John Burkardt */ { double *a_cr; double *a1; double *a2; double *a3; double *b; int i; int j; int n = 10; double *x; printf ( "\n" ); printf ( "R83V_CR_FA_TEST:\n" ); printf ( " R83V_CR_FA factors an R83V matrix using cyclic reduction;\n" ); printf ( "\n" ); printf ( " Matrix order N = %d\n", n ); printf ( " The matrix is NOT symmetric.\n" ); /* Set the matrix values. */ a1 = ( double * ) malloc ( ( n - 1 ) * sizeof ( double ) ); a2 = ( double * ) malloc ( n * sizeof ( double ) ); a3 = ( double * ) malloc ( ( n - 1 ) * sizeof ( double ) ); for ( j = 0; j < n; j++ ) { for ( i = i4_max ( 0, j - 1 ); i <= i4_min ( n - 1, j + 1 ); i++ ) { if ( j == i - 1 ) { a1[j] = ( double ) ( j + 1 ); } else if ( j == i ) { a2[j] = ( double ) ( 4 * ( j + 1 ) ); } else if ( j == i + 1 ) { a3[j-1] = ( double ) ( j + 1 ); } } } r83v_print ( n, n, a1, a2, a3, " The matrix:" ); /* Set the desired solution. */ x = r8vec_indicator1_new ( n ); /* Compute the corresponding right hand side. */ b = r83v_mv ( n, n, a1, a2, a3, x ); /* Factor the matrix. */ a_cr = r83v_cr_fa ( n, a1, a2, a3 ); /* Solve the linear system. */ x = r83v_cr_sl ( n, a_cr, b ); r8vec_print ( n, x, " Solution:" ); free ( a_cr ); free ( a1 ); free ( a2 ); free ( a3 ); free ( b ); free ( x ); return; } /******************************************************************************/ void r83v_cr_sl_test ( ) /******************************************************************************/ /* Purpose: R83V_CR_SL_TEST tests R83V_CR_SL. Licensing: This code is distributed under the MIT license. Modified: 15 February 2016 Author: John Burkardt */ { double *a_cr; double *a1; double *a2; double *a3; double *b; int i; int j; int n = 10; double *x; printf ( "\n" ); printf ( "R83V_CR_SL_TEST:\n" ); printf ( " R83V_CR_SL solves a linear system factored by R83V_CR_FA.\n" ); printf ( "\n" ); printf ( " Matrix order N = %d\n", n ); printf ( " The matrix is NOT symmetric.\n" ); /* Set the matrix values. */ a1 = ( double * ) malloc ( ( n - 1 ) * sizeof ( double ) ); a2 = ( double * ) malloc ( n * sizeof ( double ) ); a3 = ( double * ) malloc ( ( n - 1 ) * sizeof ( double ) ); for ( j = 0; j < n; j++ ) { for ( i = i4_max ( 0, j - 1 ); i <= i4_min ( n - 1, j + 1 ); i++ ) { if ( j == i - 1 ) { a1[j] = ( double ) ( j + 1 ); } else if ( j == i ) { a2[j] = ( double ) ( 4 * ( j + 1 ) ); } else if ( j == i + 1 ) { a3[j-1] = ( double ) ( j + 1 ); } } } r83v_print ( n, n, a1, a2, a3, " The matrix:" ); /* Set the desired solution. */ x = r8vec_indicator1_new ( n ); /* Compute the corresponding right hand side. */ b = r83v_mv ( n, n, a1, a2, a3, x ); /* Factor the matrix. */ a_cr = r83v_cr_fa ( n, a1, a2, a3 ); /* Solve the linear system. */ x = r83v_cr_sl ( n, a_cr, b ); r8vec_print ( n, x, " Solution:" ); free ( a_cr ); free ( a1 ); free ( a2 ); free ( a3 ); free ( b ); free ( x ); return; } /******************************************************************************/ void r83v_cr_sls_test ( ) /******************************************************************************/ /* Purpose: R83V_CR_SLS_TEST tests R83V_CR_SLS. Licensing: This code is distributed under the MIT license. Modified: 15 February 2016 Author: John Burkardt */ { double *a_cr; double *a1; double *a2; double *a3; double *b; int i; int j; int n = 5; int nb = 2; double *x; printf ( "\n" ); printf ( "R83V_CR_SLS_TEST\n" ); printf ( " R83V_CR_SLS solves multiple linear systems A*x=b1:bn with R83V matrix\n" ); printf ( " using cyclic reduction, after factorization by R83V_CR_FA.\n" ); printf ( "\n" ); printf ( " Matrix order N = %d\n", n ); printf ( " Number of linear systems = %d\n", nb ); printf ( " Demonstrate multiple system solution method.\n" ); /* Set the matrix values. */ /* Set the matrix values. */ a1 = ( double * ) malloc ( ( n - 1 ) * sizeof ( double ) ); a2 = ( double * ) malloc ( n * sizeof ( double ) ); a3 = ( double * ) malloc ( ( n - 1 ) * sizeof ( double ) ); r83v_dif2 ( n, n, a1, a2, a3 ); r83v_print ( n, n, a1, a2, a3, " System matrix:" ); /* Factor the matrix once. */ a_cr = r83v_cr_fa ( n, a1, a2, a3 ); /* Set up the linear systems. */ b = ( double * ) malloc ( n * nb * sizeof ( double ) ); for ( j = 0; j < nb; j++ ) { for ( i = 0; i < n; i++ ) { b[i+j*n] = 0.0; } } j = 0; b[n-1+j*n] = ( double ) ( n + 1 ); j = 1; b[0 +j*n] = 1.0; b[n-1+j*n] = 1.0; r8ge_print ( n, nb, b, " RHS:" ); /* Solve the linear systems. */ x = r83v_cr_sls ( n, a_cr, nb, b ); r8ge_print ( n, nb, x, " Solutions:" ); free ( a_cr ); free ( a1 ); free ( a2 ); free ( a3 ); free ( b ); free ( x ); return; } /******************************************************************************/ void r83v_dif2_test ( ) /******************************************************************************/ /* Purpose: R83V_DIF2_TEST tests R83V_DIF2. Licensing: This code is distributed under the MIT license. Modified: 13 February 2016 Author: John Burkardt */ { double *a; int ahi; double *b; int bhi; double *c; int chi; int i; int m; int n; printf ( "\n" ); printf ( "R83V_DIF2_TEST\n" ); printf ( " R83V_DIF2 sets up an R83V second difference matrix.\n" ); printf ( " We check three cases, MN.\n" ); for ( i = 1; i <= 3; i++ ) { if ( i == 1 ) { m = 3; n = 5; } else if ( i == 2 ) { m = 5; n = 5; } else if ( i == 3 ) { m = 5; n = 3; } ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a = ( double * ) malloc ( ahi * sizeof ( double ) ); b = ( double * ) malloc ( bhi * sizeof ( double ) ); c = ( double * ) malloc ( chi * sizeof ( double ) ); r83v_dif2 ( m, n, a, b, c ); r83v_print ( m, n, a, b, c, " The R83V DIF2 matrix:" ); free ( a ); free ( b ); free ( c ); } return; } /******************************************************************************/ void r83v_fs_test ( ) /******************************************************************************/ /* Purpose: R83V_FS_TEST tests R83V_FS. Licensing: This code is distributed under the MIT license. Modified: 16 February 2016 Author: John Burkardt */ { double *a1; double *a2; double *a3; double *b; int n = 10; double *x1; double *x2; printf ( "\n" ); printf ( "R83V_FS_TEST\n" ); printf ( " R83V_FS factors and solves a linear system\n" ); printf ( " for an R83V matrix.\n" ); printf ( "\n" ); printf ( " Matrix order N = %d\n", n ); /* Set the matrix values. */ a1 = ( double * ) malloc ( ( n - 1 ) * sizeof ( double ) ); a2 = ( double * ) malloc ( n * sizeof ( double ) ); a3 = ( double * ) malloc ( ( n - 1 ) * sizeof ( double ) ); r83v_dif2 ( n, n, a1, a2, a3 ); /* Set the desired solution. */ x1 = r8vec_indicator1_new ( n ); /* Compute the corresponding right hand side. */ b = r83v_mv ( n, n, a1, a2, a3, x1 ); r8vec_print ( n, b, " The right hand side:" ); /* Solve the linear system. */ x2 = r83v_fs ( n, a1, a2, a3, b ); if ( ! x2 ) { printf ( "\n" ); printf ( " R83V_FS failed.\n" ); } else { r8vec_print ( n, x2, " Solution:" ); } /* Free memory. */ free ( a1 ); free ( a2 ); free ( a3 ); free ( b ); free ( x1 ); free ( x2 ); return; } /******************************************************************************/ void r83v_gs_sl_test ( ) /******************************************************************************/ /* Purpose: R83V_GS_SL_TEST tests R83V_GS_SL. Licensing: This code is distributed under the MIT license. Modified: 15 February 2016 Author: John Burkardt */ { double *a1; double *a2; double *a3; double *b; int i; int maxit = 25; int n = 10; double *x; printf ( "\n" ); printf ( "R83V_GS_SL_TEST\n" ); printf ( " R83V_GS_SL solves a linear system using Gauss-Seidel\n" ); printf ( " iteration for an R83V matrix.\n" ); printf ( "\n" ); printf ( " Matrix order N = %d\n", n ); printf ( " Iterations per call = %d\n", maxit ); /* Set the matrix values. */ a1 = ( double * ) malloc ( ( n - 1 ) * sizeof ( double ) ); a2 = ( double * ) malloc ( n * sizeof ( double ) ); a3 = ( double * ) malloc ( ( n - 1 ) * sizeof ( double ) ); r83v_dif2 ( n, n, a1, a2, a3 ); /* Set the desired solution. */ x = r8vec_indicator1_new ( n ); /* Compute the corresponding right hand side. */ b = r83v_mv ( n, n, a1, a2, a3, x ); r8vec_print ( n, b, " The right hand side:" ); /* Set the starting solution. */ for ( i = 0; i < n; i++ ) { x[i] = 0.0; } /* Solve the linear system. */ for ( i = 1; i <= 3; i++ ) { r83v_gs_sl ( n, a1, a2, a3, b, x, maxit ); r8vec_print ( n, x, " Current estimated solution:" ); } free ( a1 ); free ( a2 ); free ( a3 ); free ( b ); free ( x ); return; } /******************************************************************************/ void r83v_indicator_test ( ) /******************************************************************************/ /* Purpose: R83V_INDICATOR_TEST tests R83V_INDICATOR. Licensing: This code is distributed under the MIT license. Modified: 13 February 2016 Author: John Burkardt */ { double *a; int ahi; double *b; int bhi; double *c; int chi; int i; int m; int n; printf ( "\n" ); printf ( "R83V_INDICATOR_TEST\n" ); printf ( " R83V_INDICATOR sets up an R83V indicator matrix.\n" ); printf ( " We check three cases, MN.\n" ); for ( i = 1; i <= 3; i++ ) { if ( i == 1 ) { m = 3; n = 5; } else if ( i == 2 ) { m = 5; n = 5; } else if ( i == 3 ) { m = 5; n = 3; } ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a = ( double * ) malloc ( ahi * sizeof ( double ) ); b = ( double * ) malloc ( bhi * sizeof ( double ) ); c = ( double * ) malloc ( chi * sizeof ( double ) ); r83v_indicator ( m, n, a, b, c ); r83v_print ( m, n, a, b, c, " The R83V indicator matrix:" ); free ( a ); free ( b ); free ( c ); } return; } /******************************************************************************/ void r83v_jac_sl_test ( ) /******************************************************************************/ /* Purpose: R83V_JAC_SL_TEST tests R83V_JAC_SL. Licensing: This code is distributed under the MIT license. Modified: 14 February 2016 Author: John Burkardt */ { double *a1; double *a2; double *a3; double *b; int i; int maxit = 25; int n = 10; double *x; printf ( "\n" ); printf ( "R83V_JAC_SL_TEST\n" ); printf ( " R83V_JAC_SL solves a linear system using Jacobi iteration,\n" ); printf ( " for an R83V matrix.\n" ); printf ( "\n" ); printf ( " Matrix order N = %d\n", n ); printf ( " Iterations per call = %d\n", maxit ); /* Set the matrix values. */ a1 = ( double * ) malloc ( ( n - 1 ) * sizeof ( double ) ); a2 = ( double * ) malloc ( n * sizeof ( double ) ); a3 = ( double * ) malloc ( ( n - 1 ) * sizeof ( double ) ); r83v_dif2 ( n, n, a1, a2, a3 ); /* Set the desired solution. */ x = r8vec_indicator1_new ( n ); /* Compute the corresponding right hand side. */ b = r83v_mv ( n, n, a1, a2, a3, x ); r8vec_print ( n, b, " The right hand side:" ); /* Set the starting solution. */ for ( i = 0; i < n; i++ ) { x[i] = 0.0; } /* Solve the linear system. */ for ( i = 1; i <= 3; i++ ) { r83v_jac_sl ( n, a1, a2, a3, b, x, maxit ); r8vec_print ( n, x, " Current estimated solution:" ); } free ( a1 ); free ( a2 ); free ( a3 ); free ( b ); free ( x ); return; } /******************************************************************************/ void r83v_mtv_test ( ) /******************************************************************************/ /* Purpose: R83V_MTV_TEST tests R83V_MTV. Licensing: This code is distributed under the MIT license. Modified: 14 February 2016 Author: John Burkardt */ { double *a_ge; double *a1; double *a2; double *a3; int ahi; double *ax; double *ax_ge; int bhi; int chi; int i; int m; int n; int seed; double *x; printf ( "\n" ); printf ( "R83V_MTV_TEST\n" ); printf ( " R83V_MTV computes b=A'*x, where A is an R83V matrix.\n" ); printf ( " We check three cases, MN.\n" ); for ( i = 1; i <= 3; i++ ) { if ( i == 1 ) { m = 3; n = 5; } else if ( i == 2 ) { m = 5; n = 5; } else if ( i == 3 ) { m = 5; n = 3; } ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a1 = ( double * ) malloc ( ahi * sizeof ( double ) ); a2 = ( double * ) malloc ( bhi * sizeof ( double ) ); a3 = ( double * ) malloc ( chi * sizeof ( double ) ); seed = 123456789; r83v_random ( m, n, &seed, a1, a2, a3 ); x = r8vec_indicator1_new ( m ); ax = r83v_mtv ( m, n, a1, a2, a3, x ); a_ge = r83v_to_r8ge ( m, n, a1, a2, a3 ); ax_ge = r8ge_mtv ( m, n, a_ge, x ); r8vec2_print ( n, ax, ax_ge, " Product comparison:" ); free ( a_ge ); free ( a1 ); free ( a2 ); free ( a3 ); free ( ax ); free ( ax_ge ); free ( x ); } return; } /******************************************************************************/ void r83v_mv_test ( ) /******************************************************************************/ /* Purpose: R83V_MV_TEST tests R83_MV. Licensing: This code is distributed under the MIT license. Modified: 14 February 2016 Author: John Burkardt */ { double *a; double *a_ge; int ahi; double *ax; double *ax_ge; double *b; int bhi; double *c; int chi; int i; int m; int n; int seed; double *x; printf ( "\n" ); printf ( "R83V_MV_TEST\n" ); printf ( " R83V_MV computes b=A*x, where A is an R83V matrix.\n" ); printf ( " We check three cases, MN.\n" ); for ( i = 1; i <= 3; i++ ) { if ( i == 1 ) { m = 3; n = 5; } else if ( i == 2 ) { m = 5; n = 5; } else if ( i == 3 ) { m = 5; n = 3; } ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a = ( double * ) malloc ( ahi * sizeof ( double ) ); b = ( double * ) malloc ( bhi * sizeof ( double ) ); c = ( double * ) malloc ( chi * sizeof ( double ) ); seed = 123456789; r83v_random ( m, n, &seed, a, b, c ); x = r8vec_indicator1_new ( n ); ax = r83v_mv ( m, n, a, b, c, x ); a_ge = r83v_to_r8ge ( m, n, a, b, c ); ax_ge = r8ge_mv ( m, n, a_ge, x ); r8vec2_print ( m, ax, ax_ge, " Product comparison:" ); free ( a ); free ( a_ge ); free ( ax ); free ( ax_ge ); free ( b ); free ( c ); free ( x ); } return; } /******************************************************************************/ void r83v_print_test ( ) /******************************************************************************/ /* Purpose: R83V_PRINT_TEST tests R83V_PRINT. Licensing: This code is distributed under the MIT license. Modified: 13 February 2016 Author: John Burkardt */ { double *a; int ahi; double *b; int bhi; double *c; int chi; int m; int n; printf ( "\n" ); printf ( "R83V_PRINT_TEST\n" ); printf ( " R83V_PRINT prints an R83V matrix.\n" ); m = 5; n = 5; ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a = ( double * ) malloc ( ahi * sizeof ( double ) ); b = ( double * ) malloc ( bhi * sizeof ( double ) ); c = ( double * ) malloc ( chi * sizeof ( double ) ); r83v_indicator ( m, n, a, b, c ); r83v_print ( m, n, a, b, c, " The R83V matrix:" ); free ( a ); free ( b ); free ( c ); return; } /******************************************************************************/ void r83v_print_some_test ( ) /******************************************************************************/ /* Purpose: R83V_PRINT_SOME_TEST tests R83V_PRINT_SOME. Licensing: This code is distributed under the MIT license. Modified: 13 February 2016 Author: John Burkardt */ { double *a; int ahi; double *b; int bhi; double *c; int chi; int m; int n; printf ( "\n" ); printf ( "R83V_PRINT_SOME_TEST\n" ); printf ( " R83V_PRINT_SOME prints some of an R83V matrix.\n" ); m = 5; n = 5; ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a = ( double * ) malloc ( ahi * sizeof ( double ) ); b = ( double * ) malloc ( bhi * sizeof ( double ) ); c = ( double * ) malloc ( chi * sizeof ( double ) ); r83v_indicator ( m, n, a, b, c ); r83v_print_some ( m, n, a, b, c, 2, 2, 5, 4, " Rows 2-5, Cols 2-4:" ); free ( a ); free ( b ); free ( c ); return; } /******************************************************************************/ void r83v_random_test ( ) /******************************************************************************/ /* Purpose: R83V_RANDOM_TEST tests R83V_RANDOM. Licensing: This code is distributed under the MIT license. Modified: 13 February 2016 Author: John Burkardt */ { double *a; int ahi; double *b; int bhi; double *c; int chi; int i; int m; int n; int seed; printf ( "\n" ); printf ( "R83V_RANDOM_TEST\n" ); printf ( " R83V_RANDOM sets up an R83V random matrix.\n" ); printf ( " We check three cases, MN.\n" ); for ( i = 1; i <= 3; i++ ) { seed = 123456789; if ( i == 1 ) { m = 3; n = 5; } else if ( i == 2 ) { m = 5; n = 5; } else if ( i == 3 ) { m = 5; n = 3; } ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a = ( double * ) malloc ( ahi * sizeof ( double ) ); b = ( double * ) malloc ( bhi * sizeof ( double ) ); c = ( double * ) malloc ( chi * sizeof ( double ) ); r83v_random ( m, n, &seed, a, b, c ); r83v_print ( m, n, a, b, c, " The R83V random matrix:" ); free ( a ); free ( b ); free ( c ); } return; } /******************************************************************************/ void r83v_res_test ( ) /******************************************************************************/ /* Purpose: R83V_RES_TEST tests R83V_RES. Licensing: This code is distributed under the MIT license. Modified: 14 February 2016 Author: John Burkardt */ { double *a; int ahi; double *ax; double *b; int bhi; double *c; int chi; int i; int m; int n; double *r; int seed; double *x; printf ( "\n" ); printf ( "R83V_RES_TEST\n" ); printf ( " R83V_RES computes b-A*x, where A is an R83V matrix.\n" ); printf ( " We check three cases, MN.\n" ); for ( i = 1; i <= 3; i++ ) { if ( i == 1 ) { m = 3; n = 5; } else if ( i == 2 ) { m = 5; n = 5; } else if ( i == 3 ) { m = 5; n = 3; } ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a = ( double * ) malloc ( ahi * sizeof ( double ) ); b = ( double * ) malloc ( bhi * sizeof ( double ) ); c = ( double * ) malloc ( chi * sizeof ( double ) ); seed = 123456789; r83v_random ( m, n, &seed, a, b, c ); x = r8vec_indicator1_new ( n ); ax = r83v_mv ( m, n, a, b, c, x ); r = r83v_res ( m, n, a, b, c, x, ax ); r8vec_print ( m, r, " Residual A*x-b:" ); free ( a ); free ( ax ); free ( b ); free ( c ); free ( r ); free ( x ); } return; } /******************************************************************************/ void r83v_to_r8ge_test ( ) /******************************************************************************/ /* Purpose: R83V_TO_R8GE_TEST tests R83V_TO_R8GE. Licensing: This code is distributed under the MIT license. Modified: 13 February 2016 Author: John Burkardt */ { double *a; double *a1; double *a2; double *a3; int ahi; int bhi; int chi; int i; int m; int n; printf ( "\n" ); printf ( "R83V_TO_R8GE_TEST\n" ); printf ( " R83V_TO_R8GE copies an R83V matrix to an R8GE matrix.\n" ); printf ( " We check three cases, MN.\n" ); for ( i = 1; i <= 3; i++ ) { if ( i == 1 ) { m = 3; n = 5; } else if ( i == 2 ) { m = 5; n = 5; } else if ( i == 3 ) { m = 5; n = 3; } ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a1 = ( double * ) malloc ( ahi * sizeof ( double ) ); a2 = ( double * ) malloc ( bhi * sizeof ( double ) ); a3 = ( double * ) malloc ( chi * sizeof ( double ) ); r83v_indicator ( m, n, a1, a2, a3 ); r83v_print ( m, n, a1, a2, a3, " R83V matrix A:" ); a = r83v_to_r8ge ( m, n, a1, a2, a3 ); r8ge_print ( m, n, a, " R8GE version of A:" ); free ( a ); free ( a1 ); free ( a2 ); free ( a3 ); } return; } /******************************************************************************/ void r83v_to_r8vec_test ( ) /******************************************************************************/ /* Purpose: R83V_TO_R8VEC_TEST tests R83V_TO_R8VEC. Licensing: This code is distributed under the MIT license. Modified: 13 February 2016 Author: John Burkardt */ { double *a; double *a1; double *a2; double *a3; int ahi; int bhi; int chi; int i; int m; int n; printf ( "\n" ); printf ( "R83V_TO_R8VEC_TEST\n" ); printf ( " R83V_TO_R8VEC copies an R83V matrix to an R8VEC.\n" ); printf ( " We check three cases, MN.\n" ); for ( i = 1; i <= 3; i++ ) { if ( i == 1 ) { m = 3; n = 5; } else if ( i == 2 ) { m = 5; n = 5; } else if ( i == 3 ) { m = 5; n = 3; } ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a1 = ( double * ) malloc ( ahi * sizeof ( double ) ); a2 = ( double * ) malloc ( bhi * sizeof ( double ) ); a3 = ( double * ) malloc ( chi * sizeof ( double ) ); r83v_indicator ( m, n, a1, a2, a3 ); r83v_print ( m, n, a1, a2, a3, " R83V matrix A:" ); a = r83v_to_r8vec ( m, n, a1, a2, a3 ); r8vec_print ( ahi + bhi + chi, a, " Vector version of A:" ); free ( a ); free ( a1 ); free ( a2 ); free ( a3 ); } return; } /******************************************************************************/ void r83v_transpose_test ( ) /******************************************************************************/ /* Purpose: R83V_TRANSPOSE_TEST tests R83V_TRANSPOSE. Licensing: This code is distributed under the MIT license. Modified: 13 February 2016 Author: John Burkardt */ { double *a1; double *a2; double *a3; int ahi; double *b1; double *b2; double *b3; int bhi; int chi; int i; int m; int n; printf ( "\n" ); printf ( "R83V_TRANSPOSE_TEST\n" ); printf ( " R83V_TRANSPOSE makes a transposed copy of an R83V matrix.\n" ); printf ( " We check three cases, MN.\n" ); for ( i = 1; i <= 3; i++ ) { if ( i == 1 ) { m = 3; n = 5; } else if ( i == 2 ) { m = 5; n = 5; } else if ( i == 3 ) { m = 5; n = 3; } ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a1 = ( double * ) malloc ( ahi * sizeof ( double ) ); a2 = ( double * ) malloc ( bhi * sizeof ( double ) ); a3 = ( double * ) malloc ( chi * sizeof ( double ) ); r83v_indicator ( m, n, a1, a2, a3 ); r83v_print ( m, n, a1, a2, a3, " R83V matrix A:" ); b1 = ( double * ) malloc ( chi * sizeof ( double ) ); b2 = ( double * ) malloc ( bhi * sizeof ( double ) ); b3 = ( double * ) malloc ( ahi * sizeof ( double ) ); r83v_transpose ( m, n, a1, a2, a3, b1, b2, b3 ); r83v_print ( n, m, b1, b2, b3, " B = copy of A:" ); free ( a1 ); free ( a2 ); free ( a3 ); free ( b1 ); free ( b2 ); free ( b3 ); } return; } /******************************************************************************/ void r83v_zeros_test ( ) /******************************************************************************/ /* Purpose: R83V_ZEROS_TEST tests R83V_ZEROS. Licensing: This code is distributed under the MIT license. Modified: 13 February 2016 Author: John Burkardt */ { double *a; int ahi; double *b; int bhi; double *c; int chi; int i; int m; int n; printf ( "\n" ); printf ( "R83V_ZEROS_TEST\n" ); printf ( " R83V_ZEROS sets up an R83V zero matrix.\n" ); printf ( " We check three cases, MN.\n" ); for ( i = 1; i <= 3; i++ ) { if ( i == 1 ) { m = 3; n = 5; } else if ( i == 2 ) { m = 5; n = 5; } else if ( i == 3 ) { m = 5; n = 3; } ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a = ( double * ) malloc ( ahi * sizeof ( double ) ); b = ( double * ) malloc ( bhi * sizeof ( double ) ); c = ( double * ) malloc ( chi * sizeof ( double ) ); r83v_zeros ( m, n, a, b, c ); r83v_print ( m, n, a, b, c, " The R83V zero matrix:" ); free ( a ); free ( b ); free ( c ); } return; } /******************************************************************************/ void r8ge_to_r83v_test ( ) /******************************************************************************/ /* Purpose: R8GE_TO_R83V_TEST tests R8GE_TO_R83V. Licensing: This code is distributed under the MIT license. Modified: 14 February 2016 Author: John Burkardt */ { double *a; double *a1; double *a2; double *a3; int ahi; int bhi; int chi; int i; int m; int n; printf ( "\n" ); printf ( "R8GE_TO_R83V_TEST\n" ); printf ( " R8GE_TO_R83V copies an R8GE matrix to an R83V matrix.\n" ); printf ( " We check three cases, MN.\n" ); for ( i = 1; i <= 3; i++ ) { if ( i == 1 ) { m = 3; n = 5; } else if ( i == 2 ) { m = 5; n = 5; } else if ( i == 3 ) { m = 5; n = 3; } a = r8ge_indicator ( m, n ); r8ge_print ( m, n, a, " R8GE matrix A:" ); ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a1 = ( double * ) malloc ( ahi * sizeof ( double ) ); a2 = ( double * ) malloc ( bhi * sizeof ( double ) ); a3 = ( double * ) malloc ( chi * sizeof ( double ) ); r8ge_to_r83v ( m, n, a, a1, a2, a3 ); r83v_print ( m, n, a1, a2, a3, " R83V copy of (some of ) matrix A:" ); free ( a ); free ( a1 ); free ( a2 ); free ( a3 ); } return; } /******************************************************************************/ void r8vec_to_r83v_test ( ) /******************************************************************************/ /* Purpose: R8VEC_TO_R83V_TEST tests R8VEC_TO_R83V. Licensing: This code is distributed under the MIT license. Modified: 15 February 2016 Author: John Burkardt */ { double *a; double *a1; double *a2; double *a3; int ahi; int bhi; int chi; int i; int m; int n; printf ( "\n" ); printf ( "R8VEC_TO_R83V_TEST\n" ); printf ( " R8VEC_TO_R83V copies an R8VEC to an R83V matrix.\n" ); printf ( " We check three cases, MN.\n" ); for ( i = 1; i <= 3; i++ ) { if ( i == 1 ) { m = 3; n = 5; } else if ( i == 2 ) { m = 5; n = 5; } else if ( i == 3 ) { m = 5; n = 3; } ahi = i4_min ( m - 1, n ); bhi = i4_min ( m, n ); chi = i4_min ( m, n - 1 ); a = r8vec_indicator1_new ( ahi + bhi + chi ); r8vec_print ( ahi + bhi + chi, a, " R8VEC:" ); a1 = ( double * ) malloc ( ahi * sizeof ( double ) ); a2 = ( double * ) malloc ( bhi * sizeof ( double ) ); a3 = ( double * ) malloc ( chi * sizeof ( double ) ); r8vec_to_r83v ( m, n, a, a1, a2, a3 ); r83v_print ( m, n, a1, a2, a3, " R83V matrix:" ); free ( a ); free ( a1 ); free ( a2 ); free ( a3 ); } return; } /******************************************************************************/ void timestamp ( ) /******************************************************************************/ /* 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: 24 September 2003 Author: John Burkardt */ { # define TIME_SIZE 40 static char time_buffer[TIME_SIZE]; const struct tm *tm; time_t now; now = time ( NULL ); tm = localtime ( &now ); strftime ( time_buffer, TIME_SIZE, "%d %B %Y %I:%M:%S %p", tm ); fprintf ( stdout, "%s\n", time_buffer ); return; # undef TIME_SIZE }