3 * Copyright (C) 2001, 2004, 2007 Brian Gough
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 3 of the License, or (at
8 * your option) any later version.
10 * This program is distributed in the hope that it will be useful, but
11 * WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 * General Public License for more details.
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
20 /* Balance a general matrix by scaling the columns
24 * where D is a diagonal matrix
29 #include <gsl/gsl_math.h>
30 #include <gsl/gsl_vector.h>
31 #include <gsl/gsl_matrix.h>
32 #include <gsl/gsl_blas.h>
34 #include <gsl/gsl_linalg.h>
37 gsl_linalg_balance_columns (gsl_matrix * A, gsl_vector * D)
39 const size_t N = A->size2;
42 if (D->size != A->size2)
44 GSL_ERROR("length of D must match second dimension of A", GSL_EINVAL);
47 gsl_vector_set_all (D, 1.0);
49 for (j = 0; j < N; j++)
51 gsl_vector_view A_j = gsl_matrix_column (A, j);
53 double s = gsl_blas_dasum(&A_j.vector);
57 if (s == 0.0 || !gsl_finite(s))
59 gsl_vector_set (D, j, f);
63 /* FIXME: we could use frexp() here */
77 gsl_vector_set (D, j, f);
81 gsl_blas_dscal(1.0/f, &A_j.vector);