1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184
|
/* eigen/symm.c
*
* Copyright (C) 2001, 2007 Brian Gough
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 3 of the License, or (at
* your option) any later version.
*
* This program is distributed in the hope that it will be useful, but
* WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
* General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program; if not, write to the Free Software
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
*/
#include <config.h>
#include <stdlib.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_vector.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_linalg.h>
#include <gsl/gsl_eigen.h>
/* Compute eigenvalues/eigenvectors of real symmetric matrix using
reduction to tridiagonal form, followed by QR iteration with
implicit shifts.
See Golub & Van Loan, "Matrix Computations" (3rd ed), Section 8.3
*/
#include "qrstep.c"
gsl_eigen_symm_workspace *
gsl_eigen_symm_alloc (const size_t n)
{
gsl_eigen_symm_workspace *w;
if (n == 0)
{
GSL_ERROR_NULL ("matrix dimension must be positive integer",
GSL_EINVAL);
}
w = ((gsl_eigen_symm_workspace *)
malloc (sizeof (gsl_eigen_symm_workspace)));
if (w == 0)
{
GSL_ERROR_NULL ("failed to allocate space for workspace", GSL_ENOMEM);
}
w->d = (double *) malloc (n * sizeof (double));
if (w->d == 0)
{
GSL_ERROR_NULL ("failed to allocate space for diagonal", GSL_ENOMEM);
}
w->sd = (double *) malloc (n * sizeof (double));
if (w->sd == 0)
{
GSL_ERROR_NULL ("failed to allocate space for subdiagonal", GSL_ENOMEM);
}
w->size = n;
return w;
}
void
gsl_eigen_symm_free (gsl_eigen_symm_workspace * w)
{
RETURN_IF_NULL (w);
free (w->sd);
free (w->d);
free (w);
}
int
gsl_eigen_symm (gsl_matrix * A, gsl_vector * eval,
gsl_eigen_symm_workspace * w)
{
if (A->size1 != A->size2)
{
GSL_ERROR ("matrix must be square to compute eigenvalues", GSL_ENOTSQR);
}
else if (eval->size != A->size1)
{
GSL_ERROR ("eigenvalue vector must match matrix size", GSL_EBADLEN);
}
else if (A->size1 != w->size)
{
GSL_ERROR ("matrix does not match workspace", GSL_EBADLEN);
}
else
{
const size_t N = A->size1;
double *const d = w->d;
double *const sd = w->sd;
size_t a, b;
/* handle special case */
if (N == 1)
{
double A00 = gsl_matrix_get (A, 0, 0);
gsl_vector_set (eval, 0, A00);
return GSL_SUCCESS;
}
/* use sd as the temporary workspace for the decomposition,
since we can discard the tau result immediately if we are not
computing eigenvectors */
{
gsl_vector_view d_vec = gsl_vector_view_array (d, N);
gsl_vector_view sd_vec = gsl_vector_view_array (sd, N - 1);
gsl_vector_view tau = gsl_vector_view_array (sd, N - 1);
gsl_linalg_symmtd_decomp (A, &tau.vector);
gsl_linalg_symmtd_unpack_T (A, &d_vec.vector, &sd_vec.vector);
}
/* Make an initial pass through the tridiagonal decomposition
to remove off-diagonal elements which are effectively zero */
chop_small_elements (N, d, sd);
/* Progressively reduce the matrix until it is diagonal */
b = N - 1;
while (b > 0)
{
if (sd[b - 1] == 0.0 || isnan(sd[b - 1]))
{
b--;
continue;
}
/* Find the largest unreduced block (a,b) starting from b
and working backwards */
a = b - 1;
while (a > 0)
{
if (sd[a - 1] == 0.0)
{
break;
}
a--;
}
{
const size_t n_block = b - a + 1;
double *d_block = d + a;
double *sd_block = sd + a;
/* apply QR reduction with implicit deflation to the
unreduced block */
qrstep (n_block, d_block, sd_block, NULL, NULL);
/* remove any small off-diagonal elements */
chop_small_elements (n_block, d_block, sd_block);
}
}
{
gsl_vector_view d_vec = gsl_vector_view_array (d, N);
gsl_vector_memcpy (eval, &d_vec.vector);
}
return GSL_SUCCESS;
}
}
|