[Top][All Lists]
[Date Prev][Date Next][Thread Prev][Thread Next][Date Index][Thread Index]
[Bug-gsl] [bug #28767] gsl_linalg_SV_decomp gives NaN when tb, tab, dt =
From: |
Brian Gough |
Subject: |
[Bug-gsl] [bug #28767] gsl_linalg_SV_decomp gives NaN when tb, tab, dt = 0 |
Date: |
Fri, 29 Jan 2010 22:05:28 +0000 |
User-agent: |
Mozilla/5.0 (X11; U; Linux i686; en-GB; rv:1.9.0.11) Gecko/2009061118 Fedora/3.0.11-1.fc9 Firefox/3.0.11 |
Follow-up Comment #1, bug #28767 (project gsl):
GSL version number: 1.9
hardware and operating system, OS details:
Linux phit1 2.6.25.20-0.5-default #1 SMP 2009-08-14 01:48:11 +0200 x86_64
x86_64 x86_64 GNU/Linux
gcc version:
gcc version 4.3.1 20080507 (prerelease) [gcc-4_3-branch revision 135036]
(SUSE Linux)
genre: linalg, SVD, gsl_linalg_SV_decomp
description of the bug behavior: NaN results in gsl_linalg_SV_decomp
short program which exercises the bug: included in this E-Mail
--------------------------------------------------------------
Dear GSL developers,
I managed to reproduce the data with a randomly generated matrix of
appropriate statistics: In this 462x421 matrix, about 90% of the rows and
43%
of the columns are completely zero. Of the other elements, 60% are 1 and the
rest is zero. Here's the test program and the output. Most of the time the
result contains NaNs; sometimes GSL's SVD algorithm fails to converge:
----- svdbug.c --------------------------------------------------
/*
GSL version: 1.9
Compile and run:
gcc `gsl-config --cflags` svdbug.c `gsl-config --libs` -o svdbug &&
./svdbug
*/
#include <stdio.h>
#include <stdlib.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_vector.h>
#define rand_double() ((double)(rand()) / RAND_MAX)
int main(int argc, char** argv) {
int m = 462, n = 421;
gsl_matrix* A = gsl_matrix_alloc(m, n);
srand(time(NULL));
int i, j;
int row_completely_zero[462];
int column_completely_zero[421];
for (i = 0; i < m; i++) {
row_completely_zero[i] = rand_double() < 0.9;
}
for (i = 0; i < n; i++) {
column_completely_zero[i] = rand_double() < 0.43;
}
int n_ones = 0, n_zeros = 0;
for (i = 0; i < m; i++) {
for (j = 0; j < n; j++) {
double a = row_completely_zero[i] == 1 || column_completely_zero[j] ==
1
? 0 : (rand_double() < 0.6 ? 1 : 0);
gsl_matrix_set(A, i, j, a);
if (a == 0) { n_zeros++; }
else if (a == 1) { n_ones++; }
}
}
printf("m = %d, n = %d, %d elements, %d zeros, %d ones, rest: %dn", m, n,
m * n, n_zeros, n_ones, m * n - n_zeros - n_ones);
gsl_matrix* V = gsl_matrix_alloc(n, n);
gsl_vector* S = gsl_vector_alloc(n);
gsl_vector* work = gsl_vector_alloc(n);
gsl_linalg_SV_decomp(A, V, S, work);
printf("first 10 singular values:n");
for (i = 0; i < 10; i++) {
if (i > 0) { printf(", "); }
printf("%g", gsl_vector_get(S, i));
}
return 0;
}
---- output ----------------------------------------------------
(phit1 ~/pj) gcc `gsl-config --cflags` svdbug.c `gsl-config --libs` -
o svdbug && ./svdbug
m = 462, n = 421, 194502 elements, 189295 zeros, 5207 ones, rest: 0
first 10 singular values:
0, 0, nan, nan, nan, nan, nan, nan, nan, nan
(phit1 ~/pj) gcc `gsl-config --cflags` svdbug.c `gsl-config --libs` -
o svdbug && ./svdbug
m = 462, n = 421, 194502 elements, 188300 zeros, 6202 ones, rest: 0
first 10 singular values:
0, nan, nan, nan, nan, nan, nan, nan, nan, nan
(phit1 ~/pj) gcc `gsl-config --cflags` svdbug.c `gsl-config --libs` -
o svdbug && ./svdbug
m = 462, n = 421, 194502 elements, 189691 zeros, 4811 ones, rest: 0
first 10 singular values:
nan, nan, nan, nan, nan, nan, nan, nan, nan, nan
(phit1 ~/pj) gcc `gsl-config --cflags` svdbug.c `gsl-config --libs` -
o svdbug && ./svdbug
m = 462, n = 421, 194502 elements, 188698 zeros, 5804 ones, rest: 0
gsl: svd.c:149: ERROR: SVD decomposition failed to converge
Default GSL error handler invoked.
[1] 29225 abort ./svdbug
(phit1 ~/pj)
Kind regards
Bruno Daniel
_______________________________________________________
Reply to this item at:
<http://savannah.gnu.org/bugs/?28767>
_______________________________________________
Message sent via/by Savannah
http://savannah.gnu.org/