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
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
|
/* eigen/gsl_eigen.h
*
* Copyright (C) 1996, 1997, 1998, 1999, 2000, 2006, 2007 Gerard Jungman, Brian Gough, Patrick Alken
*
* 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.
*/
#ifndef __GSL_EIGEN_H__
#define __GSL_EIGEN_H__
#include <gsl/gsl_vector.h>
#include <gsl/gsl_matrix.h>
#undef __BEGIN_DECLS
#undef __END_DECLS
#ifdef __cplusplus
# define __BEGIN_DECLS extern "C" {
# define __END_DECLS }
#else
# define __BEGIN_DECLS /* empty */
# define __END_DECLS /* empty */
#endif
__BEGIN_DECLS
typedef struct {
size_t size;
double * d;
double * sd;
} gsl_eigen_symm_workspace;
gsl_eigen_symm_workspace * gsl_eigen_symm_alloc (const size_t n);
void gsl_eigen_symm_free (gsl_eigen_symm_workspace * w);
int gsl_eigen_symm (gsl_matrix * A, gsl_vector * eval, gsl_eigen_symm_workspace * w);
typedef struct {
size_t size;
double * d;
double * sd;
double * gc;
double * gs;
} gsl_eigen_symmv_workspace;
gsl_eigen_symmv_workspace * gsl_eigen_symmv_alloc (const size_t n);
void gsl_eigen_symmv_free (gsl_eigen_symmv_workspace * w);
int gsl_eigen_symmv (gsl_matrix * A, gsl_vector * eval, gsl_matrix * evec, gsl_eigen_symmv_workspace * w);
typedef struct {
size_t size;
double * d;
double * sd;
double * tau;
} gsl_eigen_herm_workspace;
gsl_eigen_herm_workspace * gsl_eigen_herm_alloc (const size_t n);
void gsl_eigen_herm_free (gsl_eigen_herm_workspace * w);
int gsl_eigen_herm (gsl_matrix_complex * A, gsl_vector * eval,
gsl_eigen_herm_workspace * w);
typedef struct {
size_t size;
double * d;
double * sd;
double * tau;
double * gc;
double * gs;
} gsl_eigen_hermv_workspace;
gsl_eigen_hermv_workspace * gsl_eigen_hermv_alloc (const size_t n);
void gsl_eigen_hermv_free (gsl_eigen_hermv_workspace * w);
int gsl_eigen_hermv (gsl_matrix_complex * A, gsl_vector * eval,
gsl_matrix_complex * evec,
gsl_eigen_hermv_workspace * w);
typedef struct {
size_t size; /* matrix size */
size_t max_iterations; /* max iterations since last eigenvalue found */
size_t n_iter; /* number of iterations since last eigenvalue found */
size_t n_evals; /* number of eigenvalues found so far */
int compute_t; /* compute Schur form T = Z^t A Z */
gsl_matrix *H; /* pointer to Hessenberg matrix */
gsl_matrix *Z; /* pointer to Schur vector matrix */
} gsl_eigen_francis_workspace;
gsl_eigen_francis_workspace * gsl_eigen_francis_alloc (void);
void gsl_eigen_francis_free (gsl_eigen_francis_workspace * w);
void gsl_eigen_francis_T (const int compute_t,
gsl_eigen_francis_workspace * w);
int gsl_eigen_francis (gsl_matrix * H, gsl_vector_complex * eval,
gsl_eigen_francis_workspace * w);
int gsl_eigen_francis_Z (gsl_matrix * H, gsl_vector_complex * eval,
gsl_matrix * Z,
gsl_eigen_francis_workspace * w);
typedef struct {
size_t size; /* size of matrices */
gsl_vector *diag; /* diagonal matrix elements from balancing */
gsl_vector *tau; /* Householder coefficients */
gsl_matrix *Z; /* pointer to Z matrix */
int do_balance; /* perform balancing transformation? */
size_t n_evals; /* number of eigenvalues found */
gsl_eigen_francis_workspace *francis_workspace_p;
} gsl_eigen_nonsymm_workspace;
gsl_eigen_nonsymm_workspace * gsl_eigen_nonsymm_alloc (const size_t n);
void gsl_eigen_nonsymm_free (gsl_eigen_nonsymm_workspace * w);
void gsl_eigen_nonsymm_params (const int compute_t, const int balance,
gsl_eigen_nonsymm_workspace *w);
int gsl_eigen_nonsymm (gsl_matrix * A, gsl_vector_complex * eval,
gsl_eigen_nonsymm_workspace * w);
int gsl_eigen_nonsymm_Z (gsl_matrix * A, gsl_vector_complex * eval,
gsl_matrix * Z, gsl_eigen_nonsymm_workspace * w);
typedef struct {
size_t size; /* size of matrices */
gsl_vector *work; /* scratch workspace */
gsl_vector *work2; /* scratch workspace */
gsl_vector *work3; /* scratch workspace */
gsl_matrix *Z; /* pointer to Schur vectors */
gsl_eigen_nonsymm_workspace *nonsymm_workspace_p;
} gsl_eigen_nonsymmv_workspace;
gsl_eigen_nonsymmv_workspace * gsl_eigen_nonsymmv_alloc (const size_t n);
void gsl_eigen_nonsymmv_free (gsl_eigen_nonsymmv_workspace * w);
void gsl_eigen_nonsymmv_params (const int balance,
gsl_eigen_nonsymmv_workspace *w);
int gsl_eigen_nonsymmv (gsl_matrix * A, gsl_vector_complex * eval,
gsl_matrix_complex * evec,
gsl_eigen_nonsymmv_workspace * w);
int gsl_eigen_nonsymmv_Z (gsl_matrix * A, gsl_vector_complex * eval,
gsl_matrix_complex * evec, gsl_matrix * Z,
gsl_eigen_nonsymmv_workspace * w);
typedef struct {
size_t size; /* size of matrices */
gsl_eigen_symm_workspace *symm_workspace_p;
} gsl_eigen_gensymm_workspace;
gsl_eigen_gensymm_workspace * gsl_eigen_gensymm_alloc (const size_t n);
void gsl_eigen_gensymm_free (gsl_eigen_gensymm_workspace * w);
int gsl_eigen_gensymm (gsl_matrix * A, gsl_matrix * B,
gsl_vector * eval, gsl_eigen_gensymm_workspace * w);
int gsl_eigen_gensymm_standardize (gsl_matrix * A, const gsl_matrix * B);
typedef struct {
size_t size; /* size of matrices */
gsl_eigen_symmv_workspace *symmv_workspace_p;
} gsl_eigen_gensymmv_workspace;
gsl_eigen_gensymmv_workspace * gsl_eigen_gensymmv_alloc (const size_t n);
void gsl_eigen_gensymmv_free (gsl_eigen_gensymmv_workspace * w);
int gsl_eigen_gensymmv (gsl_matrix * A, gsl_matrix * B,
gsl_vector * eval, gsl_matrix * evec,
gsl_eigen_gensymmv_workspace * w);
typedef struct {
size_t size; /* size of matrices */
gsl_eigen_herm_workspace *herm_workspace_p;
} gsl_eigen_genherm_workspace;
gsl_eigen_genherm_workspace * gsl_eigen_genherm_alloc (const size_t n);
void gsl_eigen_genherm_free (gsl_eigen_genherm_workspace * w);
int gsl_eigen_genherm (gsl_matrix_complex * A, gsl_matrix_complex * B,
gsl_vector * eval, gsl_eigen_genherm_workspace * w);
int gsl_eigen_genherm_standardize (gsl_matrix_complex * A,
const gsl_matrix_complex * B);
typedef struct {
size_t size; /* size of matrices */
gsl_eigen_hermv_workspace *hermv_workspace_p;
} gsl_eigen_genhermv_workspace;
gsl_eigen_genhermv_workspace * gsl_eigen_genhermv_alloc (const size_t n);
void gsl_eigen_genhermv_free (gsl_eigen_genhermv_workspace * w);
int gsl_eigen_genhermv (gsl_matrix_complex * A, gsl_matrix_complex * B,
gsl_vector * eval, gsl_matrix_complex * evec,
gsl_eigen_genhermv_workspace * w);
typedef struct {
size_t size; /* size of matrices */
gsl_vector *work; /* scratch workspace */
size_t n_evals; /* number of eigenvalues found */
size_t max_iterations; /* maximum QZ iterations allowed */
size_t n_iter; /* number of iterations since last eigenvalue found */
double eshift; /* exceptional shift counter */
int needtop; /* need to compute top index? */
double atol; /* tolerance for splitting A matrix */
double btol; /* tolerance for splitting B matrix */
double ascale; /* scaling factor for shifts */
double bscale; /* scaling factor for shifts */
gsl_matrix *H; /* pointer to hessenberg matrix */
gsl_matrix *R; /* pointer to upper triangular matrix */
int compute_s; /* compute generalized Schur form S */
int compute_t; /* compute generalized Schur form T */
gsl_matrix *Q; /* pointer to left Schur vectors */
gsl_matrix *Z; /* pointer to right Schur vectors */
} gsl_eigen_gen_workspace;
gsl_eigen_gen_workspace * gsl_eigen_gen_alloc (const size_t n);
void gsl_eigen_gen_free (gsl_eigen_gen_workspace * w);
void gsl_eigen_gen_params (const int compute_s, const int compute_t,
const int balance, gsl_eigen_gen_workspace * w);
int gsl_eigen_gen (gsl_matrix * A, gsl_matrix * B,
gsl_vector_complex * alpha, gsl_vector * beta,
gsl_eigen_gen_workspace * w);
int gsl_eigen_gen_QZ (gsl_matrix * A, gsl_matrix * B,
gsl_vector_complex * alpha, gsl_vector * beta,
gsl_matrix * Q, gsl_matrix * Z,
gsl_eigen_gen_workspace * w);
typedef struct {
size_t size; /* size of matrices */
gsl_vector *work1; /* 1-norm of columns of A */
gsl_vector *work2; /* 1-norm of columns of B */
gsl_vector *work3; /* real part of eigenvector */
gsl_vector *work4; /* imag part of eigenvector */
gsl_vector *work5; /* real part of back-transformed eigenvector */
gsl_vector *work6; /* imag part of back-transformed eigenvector */
gsl_matrix *Q; /* pointer to left Schur vectors */
gsl_matrix *Z; /* pointer to right Schur vectors */
gsl_eigen_gen_workspace *gen_workspace_p;
} gsl_eigen_genv_workspace;
gsl_eigen_genv_workspace * gsl_eigen_genv_alloc (const size_t n);
void gsl_eigen_genv_free (gsl_eigen_genv_workspace * w);
int gsl_eigen_genv (gsl_matrix * A, gsl_matrix * B,
gsl_vector_complex * alpha, gsl_vector * beta,
gsl_matrix_complex * evec,
gsl_eigen_genv_workspace * w);
int gsl_eigen_genv_QZ (gsl_matrix * A, gsl_matrix * B,
gsl_vector_complex * alpha, gsl_vector * beta,
gsl_matrix_complex * evec,
gsl_matrix * Q, gsl_matrix * Z,
gsl_eigen_genv_workspace * w);
typedef enum {
GSL_EIGEN_SORT_VAL_ASC,
GSL_EIGEN_SORT_VAL_DESC,
GSL_EIGEN_SORT_ABS_ASC,
GSL_EIGEN_SORT_ABS_DESC
}
gsl_eigen_sort_t;
/* Sort eigensystem results based on eigenvalues.
* Sorts in order of increasing value or increasing
* absolute value.
*
* exceptions: GSL_EBADLEN
*/
int gsl_eigen_symmv_sort(gsl_vector * eval, gsl_matrix * evec,
gsl_eigen_sort_t sort_type);
int gsl_eigen_hermv_sort(gsl_vector * eval, gsl_matrix_complex * evec,
gsl_eigen_sort_t sort_type);
int gsl_eigen_nonsymmv_sort(gsl_vector_complex * eval,
gsl_matrix_complex * evec,
gsl_eigen_sort_t sort_type);
int gsl_eigen_gensymmv_sort (gsl_vector * eval, gsl_matrix * evec,
gsl_eigen_sort_t sort_type);
int gsl_eigen_genhermv_sort (gsl_vector * eval, gsl_matrix_complex * evec,
gsl_eigen_sort_t sort_type);
int gsl_eigen_genv_sort (gsl_vector_complex * alpha, gsl_vector * beta,
gsl_matrix_complex * evec,
gsl_eigen_sort_t sort_type);
/* Prototypes for the schur module */
int gsl_schur_gen_eigvals(const gsl_matrix *A, const gsl_matrix *B,
double *wr1, double *wr2, double *wi,
double *scale1, double *scale2);
int gsl_schur_solve_equation(double ca, const gsl_matrix *A, double z,
double d1, double d2, const gsl_vector *b,
gsl_vector *x, double *s, double *xnorm,
double smin);
int gsl_schur_solve_equation_z(double ca, const gsl_matrix *A,
gsl_complex *z, double d1, double d2,
const gsl_vector_complex *b,
gsl_vector_complex *x, double *s,
double *xnorm, double smin);
/* The following functions are obsolete: */
/* Eigensolve by Jacobi Method
*
* The data in the matrix input is destroyed.
*
* exceptions:
*/
int
gsl_eigen_jacobi(gsl_matrix * matrix,
gsl_vector * eval,
gsl_matrix * evec,
unsigned int max_rot,
unsigned int * nrot);
/* Invert by Jacobi Method
*
* exceptions:
*/
int
gsl_eigen_invert_jacobi(const gsl_matrix * matrix,
gsl_matrix * ainv,
unsigned int max_rot);
__END_DECLS
#endif /* __GSL_EIGEN_H__ */
|