123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507 |
- /*****************************************************************************
- *
- * MODULE: Grass PDE Numerical Library
- * AUTHOR(S): Soeren Gebbert, Berlin (GER) Dec 2006
- * soerengebbert <at> gmx <dot> de
- *
- * PURPOSE: Unit tests for les solving
- *
- * COPYRIGHT: (C) 2000 by the GRASS Development Team
- *
- * This program is free software under the GNU General Public
- * License (>=v2). Read the file COPYING that comes with GRASS
- * for details.
- *
- *****************************************************************************/
- #include <stdio.h>
- #include <stdlib.h>
- #include <string.h>
- #include <grass/glocale.h>
- #include <grass/gmath.h>
- #include "test_gmath_lib.h"
- #define EPSILON_DIRECT 1.0E-10
- #define EPSILON_ITER 1.0E-4
- /* prototypes */
- static int test_solvers(void);
- /* ************************************************************************* */
- /* Performe the solver unit tests ****************************************** */
- /* ************************************************************************* */
- int unit_test_solvers(void)
- {
- int sum = 0;
- G_message(_("\n++ Running solver unit tests ++"));
- sum += test_solvers();
- if (sum > 0)
- G_warning(_("\n-- Solver unit tests failure --"));
- else
- G_message(_("\n-- Solver unit tests finished successfully --"));
- return sum;
- }
- /* *************************************************************** */
- /* Test all implemented solvers for sparse and normal matrix *** */
- /* *************************************************************** */
- int test_solvers(void)
- {
- G_math_les *les;
- G_math_les *sples;
- int sum = 0;
- double val = 0.0;
- G_message("\t * testing jacobi solver with symmetric matrix\n");
- les = create_normal_symmetric_les(TEST_NUM_ROWS);
- sples = create_sparse_symmetric_les(TEST_NUM_ROWS);
- G_math_solver_jacobi(les->A, les->x, les->b, les->rows, 250, 1, 0.1e-10);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_jacobi abs %2.20f != %i", val,
- les->rows);
- sum++;
- }
- G_math_solver_sparse_jacobi(sples->Asp, sples->x, sples->b, les->rows, 250,
- 1, 0.1e-10);
- G_math_d_asum_norm(sples->x, &val, sples->rows);
- if ((val - (double)sples->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_sparse_jacobi abs %2.20f != %i", val,
- sples->rows);
- sum++;
- }
- G_math_free_les(les);
- G_math_free_les(sples);
- G_message("\t * testing jacobi solver with unsymmetric matrix\n");
- les = create_normal_unsymmetric_les(TEST_NUM_ROWS);
- sples = create_sparse_unsymmetric_les(TEST_NUM_ROWS);
- G_math_solver_jacobi(les->A, les->x, les->b, les->rows, 250, 1, 0.1e-10);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_jacobi abs %2.20f != %i", val,
- les->rows);
- sum++;
- }
- G_math_solver_sparse_jacobi(sples->Asp, sples->x, sples->b, les->rows, 250,
- 1, 0.1e-10);
- G_math_d_asum_norm(sples->x, &val, sples->rows);
- if ((val - (double)sples->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_sparse_jacobi abs %2.20f != %i", val,
- sples->rows);
- sum++;
- }
- G_math_free_les(les);
- G_math_free_les(sples);
- G_message("\t * testing gauss seidel solver with symmetric matrix\n");
- les = create_normal_symmetric_les(TEST_NUM_ROWS);
- sples = create_sparse_symmetric_les(TEST_NUM_ROWS);
- G_math_solver_gs(les->A, les->x, les->b, les->rows, 150, 1, 0.1e-9);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_gs abs %2.20f != %i", val, les->rows);
- sum++;
- }
- G_math_solver_sparse_gs(sples->Asp, sples->x, sples->b, les->rows, 150, 1,
- 0.1e-9);
- G_math_d_asum_norm(sples->x, &val, sples->rows);
- if ((val - (double)sples->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_sparse_gs abs %2.20f != %i", val,
- sples->rows);
- sum++;
- }
- G_math_free_les(les);
- G_math_free_les(sples);
- G_message("\t * testing gauss seidel solver with unsymmetric matrix\n");
- les = create_normal_unsymmetric_les(TEST_NUM_ROWS);
- sples = create_sparse_unsymmetric_les(TEST_NUM_ROWS);
- G_math_solver_gs(les->A, les->x, les->b, les->rows, 150, 1, 0.1e-9);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_gs abs %2.20f != %i", val, les->rows);
- sum++;
- }
- G_math_solver_sparse_gs(sples->Asp, sples->x, sples->b, les->rows, 150, 1,
- 0.1e-9);
- G_math_d_asum_norm(sples->x, &val, sples->rows);
- if ((val - (double)sples->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_sparse_gs abs %2.20f != %i", val,
- sples->rows);
- sum++;
- }
- G_math_free_les(les);
- G_math_free_les(sples);
- G_message("\t * testing pcg solver with symmetric bad conditioned matrix and preconditioner 3\n");
- les = create_normal_symmetric_pivot_les(TEST_NUM_ROWS);
- G_math_solver_pcg(les->A, les->x, les->b, les->rows, 250, 0.1e-9, 3);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_pcg abs %2.20f != %i", val, les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_free_les(les);
- G_message("\t * testing pcg solver with symmetric matrix and preconditioner 1\n");
- les = create_normal_symmetric_les(TEST_NUM_ROWS);
- sples = create_sparse_symmetric_les(TEST_NUM_ROWS);
- G_math_solver_pcg(les->A, les->x, les->b, les->rows, 250, 0.1e-9, 1);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_pcg abs %2.20f != %i", val, les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_solver_sparse_pcg(sples->Asp, sples->x, sples->b, les->rows, 250,
- 0.1e-9, 1);
- G_math_d_asum_norm(sples->x, &val, sples->rows);
- if ((val - (double)sples->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_sparse_pcg abs %2.20f != %i", val,
- sples->rows);
- sum++;
- }
- G_math_print_les(sples);
- G_math_free_les(les);
- G_math_free_les(sples);
- G_message("\t * testing pcg solver with symmetric matrix and preconditioner 2\n");
- les = create_normal_symmetric_les(TEST_NUM_ROWS);
- sples = create_sparse_symmetric_les(TEST_NUM_ROWS);
- G_math_solver_pcg(les->A, les->x, les->b, les->rows, 250, 0.1e-9, 2);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_pcg abs %2.20f != %i", val, les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_solver_sparse_pcg(sples->Asp, sples->x, sples->b, les->rows, 250,
- 0.1e-9, 2);
- G_math_d_asum_norm(sples->x, &val, sples->rows);
- if ((val - (double)sples->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_sparse_pcg abs %2.20f != %i", val,
- sples->rows);
- sum++;
- }
- G_math_print_les(sples);
- G_math_free_les(les);
- G_math_free_les(sples);
- G_message("\t * testing pcg solver with symmetric matrix and preconditioner 3\n");
- les = create_normal_symmetric_les(TEST_NUM_ROWS);
- sples = create_sparse_symmetric_les(TEST_NUM_ROWS);
- G_math_solver_pcg(les->A, les->x, les->b, les->rows, 250, 0.1e-9, 3);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_pcg abs %2.20f != %i", val, les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_solver_sparse_pcg(sples->Asp, sples->x, sples->b, les->rows, 250,
- 0.1e-9, 3);
- G_math_d_asum_norm(sples->x, &val, sples->rows);
- if ((val - (double)sples->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_sparse_pcg abs %2.20f != %i", val,
- sples->rows);
- sum++;
- }
- G_math_print_les(sples);
- G_math_free_les(les);
- G_math_free_les(sples);
- G_message("\t * testing cg solver with symmetric matrix\n");
- les = create_normal_symmetric_les(TEST_NUM_ROWS);
- sples = create_sparse_symmetric_les(TEST_NUM_ROWS);
- G_math_solver_cg(les->A, les->x, les->b, les->rows, 250, 0.1e-9);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_cg abs %2.20f != %i", val, les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_solver_sparse_cg(sples->Asp, sples->x, sples->b, les->rows, 250,
- 0.1e-9);
- G_math_d_asum_norm(sples->x, &val, sples->rows);
- if ((val - (double)sples->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_sparse_cg abs %2.20f != %i", val,
- sples->rows);
- sum++;
- }
- G_math_print_les(sples);
- G_math_free_les(les);
- G_math_free_les(sples);
- G_message("\t * testing cg solver with symmetric bad conditioned matrix\n");
- les = create_normal_symmetric_pivot_les(TEST_NUM_ROWS);
- G_math_solver_cg(les->A, les->x, les->b, les->rows, 250, 0.1e-9);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_cg abs %2.20f != %i", val, les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_free_les(les);
- G_message("\t * testing bicgstab solver with symmetric matrix\n");
- les = create_normal_symmetric_les(TEST_NUM_ROWS);
- sples = create_sparse_symmetric_les(TEST_NUM_ROWS);
- G_math_solver_bicgstab(les->A, les->x, les->b, les->rows, 250, 0.1e-9);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_bicgstab abs %2.20f != %i", val,
- les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_solver_sparse_bicgstab(sples->Asp, sples->x, sples->b, les->rows,
- 250, 0.1e-9);
- G_math_d_asum_norm(sples->x, &val, sples->rows);
- if ((val - (double)sples->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_sparse_bicgstab abs %2.20f != %i",
- val, sples->rows);
- sum++;
- }
- G_math_print_les(sples);
- G_math_free_les(les);
- G_math_free_les(sples);
- G_message("\t * testing bicgstab solver with unsymmetric matrix\n");
- les = create_normal_unsymmetric_les(TEST_NUM_ROWS);
- sples = create_sparse_unsymmetric_les(TEST_NUM_ROWS);
- G_math_solver_bicgstab(les->A, les->x, les->b, les->rows, 250, 0.1e-9);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_bicgstab abs %2.20f != %i", val,
- les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_solver_sparse_bicgstab(sples->Asp, sples->x, sples->b, les->rows,
- 250, 0.1e-9);
- G_math_d_asum_norm(sples->x, &val, sples->rows);
- if ((val - (double)les->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_sparse_bicgstab abs %2.20f != %i",
- val, sples->rows);
- sum++;
- }
- G_math_print_les(sples);
- G_math_free_les(les);
- G_math_free_les(sples);
- G_message("\t * testing gauss elimination solver with symmetric matrix\n");
- les = create_normal_symmetric_les(TEST_NUM_ROWS);
- G_math_solver_gauss(les->A, les->x, les->b, les->rows);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_DIRECT)
- {
- G_warning("Error in G_math_solver_gauss abs %2.20f != %i", val,
- les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_free_les(les);
- G_message("\t * testing lu decomposition solver with symmetric matrix\n");
- les = create_normal_symmetric_les(TEST_NUM_ROWS);
- G_math_solver_lu(les->A, les->x, les->b, les->rows);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_DIRECT)
- {
- G_warning("Error in G_math_solver_gauss abs %2.20f != %i", val,
- les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_free_les(les);
- G_message("\t * testing gauss elimination solver with unsymmetric matrix\n");
- les = create_normal_unsymmetric_les(TEST_NUM_ROWS);
- G_math_solver_gauss(les->A, les->x, les->b, les->rows);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_DIRECT)
- {
- G_warning("Error in G_math_solver_gauss abs %2.20f != %i", val,
- les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_free_les(les);
- G_message("\t * testing lu decomposition solver with unsymmetric matrix\n");
- les = create_normal_unsymmetric_les(TEST_NUM_ROWS);
- G_math_solver_lu(les->A, les->x, les->b, les->rows);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_DIRECT)
- {
- G_warning("Error in G_math_solver_gauss abs %2.20f != %i", val,
- les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_free_les(les);
- G_message("\t * testing gauss elimination solver with symmetric bad conditioned matrix\n");
- les = create_normal_symmetric_pivot_les(TEST_NUM_ROWS);
- G_math_solver_gauss(les->A, les->x, les->b, les->rows);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_DIRECT)
- {
- G_warning("Error in G_math_solver_gauss abs %2.20f != %i", val,
- les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_free_les(les);
- G_message("\t * testing lu decomposition solver with symmetric bad conditioned matrix\n");
- les = create_normal_symmetric_pivot_les(TEST_NUM_ROWS);
- G_math_solver_lu(les->A, les->x, les->b, les->rows);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_DIRECT)
- {
- G_warning("Error in G_math_solver_lu abs %2.20f != %i", val,
- les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_free_les(les);
- G_message("\t * testing cholesky decomposition solver with symmetric matrix\n");
- les = create_normal_symmetric_les(TEST_NUM_ROWS);
- /*cholesky*/G_math_solver_cholesky(les->A, les->x, les->b, les->rows,
- les->rows);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_DIRECT)
- {
- G_warning("Error in G_math_solver_solver_cholesky abs %2.20f != %i", val,
- les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_free_les(les);
- G_message("\t * testing cholesky band decomposition solver with symmetric band matrix 1\n");
- les = create_normal_symmetric_les(TEST_NUM_ROWS);
- G_math_print_les(les);
- /* Create a band matrix*/
- G_message("\t * Creating symmetric band matrix\n");
- les->A = G_math_matrix_to_sband_matrix(les->A, les->rows, les->rows);
- G_math_print_les(les);
- /*cholesky*/G_math_solver_cholesky_sband(les->A, les->x, les->b, les->rows,les->rows);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_DIRECT)
- {
- G_warning("Error in G_math_solver_solver_cholesky_band abs %2.20f != %i", val,
- les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_free_les(les);
- G_message("\t * testing cholesky band decomposition solver with symmetric band matrix 2\n");
- les = create_symmetric_band_les(TEST_NUM_ROWS);
- G_math_print_les(les);
- /*cholesky*/G_math_solver_cholesky_sband(les->A, les->x, les->b, les->rows,les->rows);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_DIRECT)
- {
- G_warning("Error in G_math_solver_solver_cholesky_band abs %2.20f != %i", val,
- les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_free_les(les);
- G_message("\t * testing cg solver with symmetric band matrix\n");
- les = create_symmetric_band_les(TEST_NUM_ROWS);
- G_math_solver_cg_sband(les->A, les->x, les->b, les->rows, les->rows, 250, 0.1e-9);
- G_math_d_asum_norm(les->x, &val, les->rows);
- if ((val - (double)les->rows) > EPSILON_ITER)
- {
- G_warning("Error in G_math_solver_cg_sband abs %2.20f != %i", val, les->rows);
- sum++;
- }
- G_math_print_les(les);
- G_math_free_les(les);
- return sum;
- }
|