3 * Copyright (C) 2001, 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.
22 #include <gsl/gsl_math.h>
23 #include <gsl/gsl_vector.h>
24 #include <gsl/gsl_matrix.h>
25 #include <gsl/gsl_linalg.h>
26 #include <gsl/gsl_eigen.h>
28 /* Compute eigenvalues of complex hermitian matrix using reduction to
29 real symmetric tridiagonal form, followed by QR iteration with
32 See Golub & Van Loan, "Matrix Computations" (3rd ed), Section 8.3 */
36 gsl_eigen_herm_workspace *
37 gsl_eigen_herm_alloc (const size_t n)
39 gsl_eigen_herm_workspace * w ;
43 GSL_ERROR_NULL ("matrix dimension must be positive integer", GSL_EINVAL);
46 w = (gsl_eigen_herm_workspace *) malloc (sizeof(gsl_eigen_herm_workspace));
50 GSL_ERROR_NULL ("failed to allocate space for workspace", GSL_ENOMEM);
53 w->d = (double *) malloc (n * sizeof (double));
57 GSL_ERROR_NULL ("failed to allocate space for diagonal", GSL_ENOMEM);
60 w->sd = (double *) malloc (n * sizeof (double));
64 GSL_ERROR_NULL ("failed to allocate space for subdiagonal", GSL_ENOMEM);
67 w->tau = (double *) malloc (2 * n * sizeof (double));
71 GSL_ERROR_NULL ("failed to allocate space for tau", GSL_ENOMEM);
80 gsl_eigen_herm_free (gsl_eigen_herm_workspace * w)
89 gsl_eigen_herm (gsl_matrix_complex * A, gsl_vector * eval,
90 gsl_eigen_herm_workspace * w)
92 if (A->size1 != A->size2)
94 GSL_ERROR ("matrix must be square to compute eigenvalues", GSL_ENOTSQR);
96 else if (eval->size != A->size1)
98 GSL_ERROR ("eigenvalue vector must match matrix size", GSL_EBADLEN);
102 const size_t N = A->size1;
103 double *const d = w->d;
104 double *const sd = w->sd;
108 /* handle special case */
112 gsl_complex A00 = gsl_matrix_complex_get (A, 0, 0);
113 gsl_vector_set (eval, 0, GSL_REAL(A00));
118 gsl_vector_view d_vec = gsl_vector_view_array (d, N);
119 gsl_vector_view sd_vec = gsl_vector_view_array (sd, N - 1);
120 gsl_vector_complex_view tau_vec = gsl_vector_complex_view_array (w->tau, N-1);
121 gsl_linalg_hermtd_decomp (A, &tau_vec.vector);
122 gsl_linalg_hermtd_unpack_T (A, &d_vec.vector, &sd_vec.vector);
125 /* Make an initial pass through the tridiagonal decomposition
126 to remove off-diagonal elements which are effectively zero */
128 chop_small_elements (N, d, sd);
130 /* Progressively reduce the matrix until it is diagonal */
136 if (sd[b - 1] == 0.0 || isnan(sd[b - 1]))
142 /* Find the largest unreduced block (a,b) starting from b
143 and working backwards */
149 if (sd[a - 1] == 0.0)
157 const size_t n_block = b - a + 1;
158 double *d_block = d + a;
159 double *sd_block = sd + a;
161 /* apply QR reduction with implicit deflation to the
164 qrstep (n_block, d_block, sd_block, NULL, NULL);
166 /* remove any small off-diagonal elements */
168 chop_small_elements (n_block, d_block, sd_block);
173 gsl_vector_view d_vec = gsl_vector_view_array (d, N);
174 gsl_vector_memcpy (eval, &d_vec.vector);