NEURON
symmeig.c
Go to the documentation of this file.
1 #include <../../nrnconf.h>
2 
3 /**************************************************************************
4 **
5 ** Copyright (C) 1993 David E. Steward & Zbigniew Leyk, all rights reserved.
6 **
7 ** Meschach Library
8 **
9 ** This Meschach Library is provided "as is" without any express
10 ** or implied warranty of any kind with respect to this software.
11 ** In particular the authors shall not be liable for any direct,
12 ** indirect, special, incidental or consequential damages arising
13 ** in any way from use of the software.
14 **
15 ** Everyone is granted permission to copy, modify and redistribute this
16 ** Meschach Library, provided:
17 ** 1. All copies contain this copyright notice.
18 ** 2. All modified copies shall carry a notice stating who
19 ** made the last modification and the date of such modification.
20 ** 3. No charge is made for this software or works derived from it.
21 ** This clause shall not be construed as constraining other software
22 ** distributed on the same medium as this software, nor is a
23 ** distribution fee considered a charge.
24 **
25 ***************************************************************************/
26 
27 
28 /*
29  File containing routines for symmetric eigenvalue problems
30 */
31 
32 #include <stdio.h>
33 #include "matrix.h"
34 #include "matrix2.h"
35 #include <math.h>
36 
37 
38 static char rcsid[] = "symmeig.c,v 1.1 1997/12/04 17:55:57 hines Exp";
39 
40 
41 
42 #define SQRT2 1.4142135623730949
43 #define sgn(x) ( (x) >= 0 ? 1 : -1 )
44 
45 /* trieig -- finds eigenvalues of symmetric tridiagonal matrices
46  -- matrix represented by a pair of vectors a (diag entries)
47  and b (sub- & super-diag entries)
48  -- eigenvalues in a on return */
49 VEC *trieig(a,b,Q)
50 VEC *a, *b;
51 MAT *Q;
52 {
53  int i, i_min, i_max, n, split;
54  Real *a_ve, *b_ve;
55  Real b_sqr, bk, ak1, bk1, ak2, bk2, z;
56  Real c, c2, cs, s, s2, d, mu;
57 
58  if ( ! a || ! b )
59  error(E_NULL,"trieig");
60  if ( a->dim != b->dim + 1 || ( Q && Q->m != a->dim ) )
61  error(E_SIZES,"trieig");
62  if ( Q && Q->m != Q->n )
63  error(E_SQUARE,"trieig");
64 
65  n = a->dim;
66  a_ve = a->ve; b_ve = b->ve;
67 
68  i_min = 0;
69  while ( i_min < n ) /* outer while loop */
70  {
71  /* find i_max to suit;
72  submatrix i_min..i_max should be irreducible */
73  i_max = n-1;
74  for ( i = i_min; i < n-1; i++ )
75  if ( b_ve[i] == 0.0 )
76  { i_max = i; break; }
77  if ( i_max <= i_min )
78  {
79  /* printf("# i_min = %d, i_max = %d\n",i_min,i_max); */
80  i_min = i_max + 1;
81  continue; /* outer while loop */
82  }
83 
84  /* printf("# i_min = %d, i_max = %d\n",i_min,i_max); */
85 
86  /* repeatedly perform QR method until matrix splits */
87  split = FALSE;
88  while ( ! split ) /* inner while loop */
89  {
90 
91  /* find Wilkinson shift */
92  d = (a_ve[i_max-1] - a_ve[i_max])/2;
93  b_sqr = b_ve[i_max-1]*b_ve[i_max-1];
94  mu = a_ve[i_max] - b_sqr/(d + sgn(d)*sqrt(d*d+b_sqr));
95  /* printf("# Wilkinson shift = %g\n",mu); */
96 
97  /* initial Givens' rotation */
98  givens(a_ve[i_min]-mu,b_ve[i_min],&c,&s);
99  s = -s;
100  /* printf("# c = %g, s = %g\n",c,s); */
101  if ( fabs(c) < SQRT2 )
102  { c2 = c*c; s2 = 1-c2; }
103  else
104  { s2 = s*s; c2 = 1-s2; }
105  cs = c*s;
106  ak1 = c2*a_ve[i_min]+s2*a_ve[i_min+1]-2*cs*b_ve[i_min];
107  bk1 = cs*(a_ve[i_min]-a_ve[i_min+1]) +
108  (c2-s2)*b_ve[i_min];
109  ak2 = s2*a_ve[i_min]+c2*a_ve[i_min+1]+2*cs*b_ve[i_min];
110  bk2 = ( i_min < i_max-1 ) ? c*b_ve[i_min+1] : 0.0;
111  z = ( i_min < i_max-1 ) ? -s*b_ve[i_min+1] : 0.0;
112  a_ve[i_min] = ak1;
113  a_ve[i_min+1] = ak2;
114  b_ve[i_min] = bk1;
115  if ( i_min < i_max-1 )
116  b_ve[i_min+1] = bk2;
117  if ( Q )
118  rot_cols(Q,i_min,i_min+1,c,-s,Q);
119  /* printf("# z = %g\n",z); */
120  /* printf("# a [temp1] =\n"); v_output(a); */
121  /* printf("# b [temp1] =\n"); v_output(b); */
122 
123  for ( i = i_min+1; i < i_max; i++ )
124  {
125  /* get Givens' rotation for sub-block -- k == i-1 */
126  givens(b_ve[i-1],z,&c,&s);
127  s = -s;
128  /* printf("# c = %g, s = %g\n",c,s); */
129 
130  /* perform Givens' rotation on sub-block */
131  if ( fabs(c) < SQRT2 )
132  { c2 = c*c; s2 = 1-c2; }
133  else
134  { s2 = s*s; c2 = 1-s2; }
135  cs = c*s;
136  bk = c*b_ve[i-1] - s*z;
137  ak1 = c2*a_ve[i]+s2*a_ve[i+1]-2*cs*b_ve[i];
138  bk1 = cs*(a_ve[i]-a_ve[i+1]) +
139  (c2-s2)*b_ve[i];
140  ak2 = s2*a_ve[i]+c2*a_ve[i+1]+2*cs*b_ve[i];
141  bk2 = ( i+1 < i_max ) ? c*b_ve[i+1] : 0.0;
142  z = ( i+1 < i_max ) ? -s*b_ve[i+1] : 0.0;
143  a_ve[i] = ak1; a_ve[i+1] = ak2;
144  b_ve[i] = bk1;
145  if ( i < i_max-1 )
146  b_ve[i+1] = bk2;
147  if ( i > i_min )
148  b_ve[i-1] = bk;
149  if ( Q )
150  rot_cols(Q,i,i+1,c,-s,Q);
151  /* printf("# a [temp2] =\n"); v_output(a); */
152  /* printf("# b [temp2] =\n"); v_output(b); */
153  }
154 
155  /* test to see if matrix should be split */
156  for ( i = i_min; i < i_max; i++ )
157  if ( fabs(b_ve[i]) < MACHEPS*
158  (fabs(a_ve[i])+fabs(a_ve[i+1])) )
159  { b_ve[i] = 0.0; split = TRUE; }
160 
161  /* printf("# a =\n"); v_output(a); */
162  /* printf("# b =\n"); v_output(b); */
163  }
164  }
165 
166  return a;
167 }
168 
169 /* symmeig -- computes eigenvalues of a dense symmetric matrix
170  -- A **must** be symmetric on entry
171  -- eigenvalues stored in out
172  -- Q contains orthogonal matrix of eigenvectors
173  -- returns vector of eigenvalues */
174 VEC *symmeig(A,Q,out)
175 MAT *A, *Q;
176 VEC *out;
177 {
178  int i;
179  static MAT *tmp = MNULL;
180  static VEC *b = VNULL, *diag = VNULL, *beta = VNULL;
181 
182  if ( ! A )
183  error(E_NULL,"symmeig");
184  if ( A->m != A->n )
185  error(E_SQUARE,"symmeig");
186  if ( ! out || out->dim != A->m )
187  out = v_resize(out,A->m);
188 
189  tmp = m_resize(tmp,A->m,A->n);
190  tmp = m_copy(A,tmp);
191  b = v_resize(b,A->m - 1);
192  diag = v_resize(diag,(u_int)A->m);
193  beta = v_resize(beta,(u_int)A->m);
194  MEM_STAT_REG(tmp,TYPE_MAT);
197  MEM_STAT_REG(beta,TYPE_VEC);
198 
199  Hfactor(tmp,diag,beta);
200  if ( Q )
201  makeHQ(tmp,diag,beta,Q);
202 
203  for ( i = 0; i < A->m - 1; i++ )
204  {
205  out->ve[i] = tmp->me[i][i];
206  b->ve[i] = tmp->me[i][i+1];
207  }
208  out->ve[i] = tmp->me[i][i];
209  trieig(out,b,Q);
210 
211  return out;
212 }
213 
#define c
#define TRUE
Definition: err.c:57
#define FALSE
Definition: err.c:56
#define E_SQUARE
Definition: err.h:103
#define error(err_num, fn_name)
Definition: err.h:73
#define E_NULL
Definition: err.h:102
#define E_SIZES
Definition: err.h:95
#define diag(s)
Definition: fmenu.cpp:192
void givens(double x, double y, Real *c, Real *s)
Definition: givens.c:47
MAT * rot_cols(MAT *mat, u_int i, u_int k, double c, double s, MAT *out)
Definition: givens.c:114
MAT * makeHQ(MAT *H, VEC *diag, VEC *beta, MAT *Qout)
Definition: hessen.c:86
MAT * Hfactor(MAT *A, VEC *diag, VEC *beta)
Definition: hessen.c:45
static Object ** v_resize(void *v)
Definition: ivocvect.cpp:1383
#define Real
Definition: machine.h:189
uint32_t u_int
Definition: machine.h:38
#define MACHEPS
Definition: machine.h:213
static Object ** m_resize(void *v)
Definition: matrix.cpp:191
#define VNULL
Definition: matrix.h:631
#define m_copy(in, out)
Definition: matrix.h:403
#define MNULL
Definition: matrix.h:632
#define i
Definition: md1redef.h:12
#define TYPE_VEC
Definition: meminfo.h:52
#define TYPE_MAT
Definition: meminfo.h:49
#define MEM_STAT_REG(var, type)
Definition: meminfo.h:133
sqrt
Definition: extdef.h:3
fabs
Definition: extdef.h:3
#define A(i)
Definition: multisplit.cpp:63
int const size_t const size_t n
Definition: nrngsl.h:11
Definition: matrix.h:73
u_int n
Definition: matrix.h:74
u_int m
Definition: matrix.h:74
Real ** me
Definition: matrix.h:76
Definition: matrix.h:67
u_int dim
Definition: matrix.h:68
Real * ve
Definition: matrix.h:69
#define SQRT2
Definition: symmeig.c:42
VEC * trieig(VEC *a, VEC *b, MAT *Q)
Definition: symmeig.c:49
VEC * symmeig(MAT *A, MAT *Q, VEC *out)
Definition: symmeig.c:174
#define sgn(x)
Definition: symmeig.c:43
static char rcsid[]
Definition: symmeig.c:38