ADMB Documentation  -a65f1c97
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
eigen.cpp
Go to the documentation of this file.
1 
8 //#define eigen_vectors
9 
10 #include <fvar.hpp>
11 #ifndef OPT_LIB
12  #include <cassert>
13  #include <climits>
14 #endif
15 
16 void tri_dag(const dmatrix& m, const dvector& d,const dvector& e);
17 void get_eigen(const dvector& d, const dvector& e, const dmatrix& z);
18 
24 {
25  if (m.rowsize()!=m.colsize())
26  {
27  cerr << "error -- non square matrix passed to "
28  "dvector eigen(const dmatrix& m)\n";
29  ad_exit(1);
30  }
31  dmatrix m1=symmetrize(m);
32  m1.colshift(1); // set minimum column and row indices to 1
33  m1.rowshift(1);
34 #if !defined(OPT_LIB) && (__cplusplus >= 201103L)
35  int n = [](unsigned int rowsize) -> int
36  {
37  assert(rowsize <= INT_MAX);
38  return static_cast<int>(rowsize);
39  } (m1.rowsize());
40 #else
41  int n = static_cast<int>(m1.rowsize());
42 #endif
43  dvector diag(1,n);
44  dvector off_diag(1,n);
45 
46  tri_dag(m1,diag,off_diag);
47 
48  get_eigen(diag,off_diag,m1); // eigenvalues are returned in diag
49  // eigenvalues are returned in columns of z
50  return diag;
51 }
52 
63 void tri_dag(const dmatrix& _m, const dvector& _d, const dvector& _e)
64 {
65  dvector& d = (dvector&) _d;
66  dvector& e = (dvector&) _e;
67  dmatrix& m = (dmatrix&) _m;
68  if (m.rowsize() != m.colsize())
69  {
70  cerr << "Error -- non square matrix passed to "
71  "void tridag(const dmatrix& m)\n";
72  ad_exit(1);
73  }
74  if (m.rowsize() != d.size() || m.rowsize() != e.size()
75  || d.indexmin() != 1 || e.indexmin() !=1 )
76  {
77  cerr <<"Error -- incorrect vector size passed to "
78  "void tridag(const dmatrix& m)\n";
79  ad_exit(1);
80  }
81 #if !defined(OPT_LIB) && (__cplusplus >= 201103L)
82  int n = [](unsigned int rowsize) -> int
83  {
84  assert(rowsize <= INT_MAX);
85  return static_cast<int>(rowsize);
86  } (m.rowsize());
87 #else
88  int n = static_cast<int>(m.rowsize());
89 #endif
90  int l,j,i;
91  double scale,hh,h,g,f;
92 
93  for (i=n;i>=2;i--)
94  {
95  l=i-1;
96  h=scale=0.0;
97  if (l > 1)
98  {
99  for (int k=1;k<=l;k++)
100  scale += fabs(m[i][k]);
101  if (scale == 0.0)
102  e[i]=m[i][l];
103  else
104  {
105  for (int k=1;k<=l;k++)
106  {
107  m[i][k] /= scale;
108  h += m[i][k]*m[i][k];
109  }
110  f=m[i][l];
111  g = f>0 ? -sqrt(h) : sqrt(h);
112  e[i]=scale*g;
113  h -= f*g;
114  m[i][l]=f-g;
115  f=0.0;
116  for (j=1;j<=l;j++)
117  {
118  #ifdef EIGEN_VECTORS
119  /* Next statement can be omitted if eigenvectors not wanted */
120  m[j][i]=m[i][j]/h;
121  #endif
122  g=0.0;
123  for (int k=1;k<=j;k++)
124  g += m[j][k]*m[i][k];
125  for (int k=j+1;k<=l;k++)
126  g += m[k][j]*m[i][k];
127  e[j]=g/h;
128  f += e[j]*m[i][j];
129  }
130  hh=f/(h+h);
131  for (j=1;j<=l;j++)
132  {
133  f=m[i][j];
134  e[j]=g=e[j]-hh*f;
135  for (int k=1;k<=j;k++)
136  m[j][k] -= (f*e[k]+g*m[i][k]);
137  }
138  }
139  }
140  else
141  {
142  e[i]=m[i][l];
143  }
144  d[i]=h;
145  }
146  /* Next statement can be omitted if eigenvectors not wanted */
147  d[1]=0.0;
148  e[1]=0.0;
149  /* Contents of this loop can be omitted if eigenvectors not
150  wanted except for statement d[i]=a[i][i]; */
151  #ifdef EIGEN_VECTORS
152  for (i=1;i<=n;i++)
153  {
154  l=i-1;
155  if (d[i])
156  {
157  for (j=1;j<=l;j++)
158  {
159  g=0.0;
160  for (int k=1;k<=l;k++)
161  g += m[i][k]*m[k][j];
162  for (int k=1;k<=l;k++)
163  m[k][j] -= g*m[k][i];
164  }
165  }
166  d[i]=m[i][i];
167  m[i][i]=1.0;
168  for (j=1;j<=l;j++) m[j][i]=m[i][j]=0.0;
169  }
170  #else
171  for (i=1;i<=n;i++)
172  {
173  d[i]=m[i][i];
174  }
175  #endif
176 }
177 
178 double SIGN(const double x, double y)
179 {
180  if (y<0)
181  {
182  return -fabs(x);
183  }
184  else
185  {
186  return fabs(x);
187  }
188 }
189 //#define SIGN(a,b) ((b)<0 ? -fabs(a) : fabs(a))
190 
201 void get_eigen(const dvector& _d, const dvector& _e, const dmatrix& _z)
202 {
203  dvector& d = (dvector&) _d;
204  dvector& e = (dvector&) _e;
205 #ifdef EIGEN_VECTORS
206  dmatrix& z = (dmatrix&) _z;
207 #endif
208  int max_iterations=30;
209 #if !defined(OPT_LIB) && (__cplusplus >= 201103L)
210  int n = [](unsigned int size) -> int
211  {
212  assert(size <= INT_MAX);
213  return static_cast<int>(size);
214  } (d.size());
215 #else
216  int n = static_cast<int>(d.size());
217 #endif
218  max_iterations+=10*(n/100);
219  int m,l,iter,i;
220  double s,r,p,g,f,dd,c,b;
221 
222  for (i=2;i<=n;i++) e[i-1]=e[i];
223  e[n]=0.0;
224  for (l=1;l<=n;l++) {
225  iter=0;
226  do {
227  for (m=l;m<=n-1;m++) {
228  dd=fabs(d[m])+fabs(d[m+1]);
229  if (fabs(e[m])+dd == dd) break;
230  }
231  if (m != l)
232  {
233  if (iter++ == max_iterations)
234  {
235  cerr << "Maximum number of iterations exceeded in"
236  " dvector eigen(const dmatrix& m)\n";
237  ad_exit(1);
238  }
239  g=(d[l+1]-d[l])/(2.0*e[l]);
240  r=sqrt((g*g)+1.0);
241  g=d[m]-d[l]+e[l]/(g+SIGN(r,g));
242  s=c=1.0;
243  p=0.0;
244  for (i=m-1;i>=l;i--) {
245  f=s*e[i];
246  b=c*e[i];
247  if (fabs(f) >= fabs(g)) {
248  c=g/f;
249  r=sqrt((c*c)+1.0);
250  e[i+1]=f*r;
251  c *= (s=1.0/r);
252  } else {
253  s=f/g;
254  r=sqrt((s*s)+1.0);
255  e[i+1]=g*r;
256  s *= (c=1.0/r);
257  }
258  g=d[i+1]-p;
259  r=(d[i]-g)*s+2.0*c*b;
260  p=s*r;
261  d[i+1]=g+p;
262  g=c*r-b;
263  /* Next loop can be omitted if eigenvectors not wanted */
264  #ifdef EIGEN_VECTORS
265  for (int k=1;k<=n;k++)
266  {
267  f=z[k][i+1];
268  z[k][i+1]=s*z[k][i]+c*f;
269  z[k][i]=c*z[k][i]-s*f;
270  }
271  #endif
272  }
273  d[l]=d[l]-p;
274  e[l]=g;
275  e[m]=0.0;
276  }
277  } while (m != l);
278  }
279 }
280 
292 {
293  dvector& d = (dvector&) _d;
294  dvector& e = (dvector&) _e;
295 
296  int max_iterations=30;
297 #if !defined(OPT_LIB) && (__cplusplus >= 201103L)
298  int n = [](unsigned int size) -> int
299  {
300  assert(size <= INT_MAX);
301  return static_cast<int>(size);
302  } (d.size());
303 #else
304  int n = static_cast<int>(d.size());
305 #endif
306  max_iterations+=10*(n/100);
307  int m,l,iter,i;
308  double s,r,p,g,f,dd,c,b;
309 
310  for (i=2;i<=n;i++) e[i-1]=e[i];
311  e[n]=0.0;
312  for (l=1;l<=n;l++) {
313  iter=0;
314  do {
315  for (m=l;m<=n-1;m++) {
316  dd=fabs(d[m])+fabs(d[m+1]);
317  if (fabs(e[m])+dd == dd) break;
318  }
319  if (m != l)
320  {
321  if (iter++ == max_iterations)
322  {
323  cerr << "Maximum number of iterations exceeded in"
324  " dvector eigen(const dmatrix& m)\n";
325  ad_exit(1);
326  }
327  g=(d[l+1]-d[l])/(2.0*e[l]);
328  r=sqrt((g*g)+1.0);
329  g=d[m]-d[l]+e[l]/(g+SIGN(r,g));
330  s=c=1.0;
331  p=0.0;
332  for (i=m-1;i>=l;i--) {
333  f=s*e[i];
334  b=c*e[i];
335  if (fabs(f) >= fabs(g)) {
336  c=g/f;
337  r=sqrt((c*c)+1.0);
338  e[i+1]=f*r;
339  c *= (s=1.0/r);
340  } else {
341  s=f/g;
342  r=sqrt((s*s)+1.0);
343  e[i+1]=g*r;
344  s *= (c=1.0/r);
345  }
346  g=d[i+1]-p;
347  r=(d[i]-g)*s+2.0*c*b;
348  p=s*r;
349  d[i+1]=g+p;
350  g=c*r-b;
351  }
352  d[l]=d[l]-p;
353  e[l]=g;
354  e[m]=0.0;
355  }
356  } while (m != l);
357  }
358  return d;
359 }
360 
373  const dmatrix& _z)
374  //eigenvectors are returned in z
375 {
376  dvector& xd = (dvector&) _d;
377  dvector& xe = (dvector&) _e;
378  dmatrix& z = (dmatrix&) _z;
379 
380  dvector d(xd.indexmin(),xd.indexmax());
381  dvector e(xe.indexmin(),xe.indexmax());
382 
383  d=xd;
384  e=xe;
385 
386  int max_iterations=30;
387 #if !defined(OPT_LIB) && (__cplusplus >= 201103L)
388  int n = [](unsigned int size) -> int
389  {
390  assert(size <= INT_MAX);
391  return static_cast<int>(size);
392  } (d.size());
393 #else
394  int n = static_cast<int>(d.size());
395 #endif
396  max_iterations+=10*(n/100);
397  int m,l,iter,i;
398  double s,r,p,g,f,dd,c,b;
399 
400  for (i=2;i<=n;i++) e[i-1]=e[i];
401  e[n]=0.0;
402  for (l=1;l<=n;l++) {
403  iter=0;
404  do {
405  for (m=l;m<=n-1;m++) {
406  dd=fabs(d[m])+fabs(d[m+1]);
407  if (fabs(e[m])+dd == dd) break;
408  }
409  if (m != l)
410  {
411  if (iter++ == max_iterations)
412  {
413  cerr << "Maximum number of iterations exceeded in"
414  " dvector eigen(const dmatrix& m)\n";
415  ad_exit(1);
416  }
417  g=(d[l+1]-d[l])/(2.0*e[l]);
418  r=sqrt((g*g)+1.0);
419  g=d[m]-d[l]+e[l]/(g+SIGN(r,g));
420  s=c=1.0;
421  p=0.0;
422  for (i=m-1;i>=l;i--) {
423  f=s*e[i];
424  b=c*e[i];
425  if (fabs(f) >= fabs(g)) {
426  c=g/f;
427  r=sqrt((c*c)+1.0);
428  e[i+1]=f*r;
429  c *= (s=1.0/r);
430  } else {
431  s=f/g;
432  r=sqrt((s*s)+1.0);
433  e[i+1]=g*r;
434  s *= (c=1.0/r);
435  }
436  g=d[i+1]-p;
437  r=(d[i]-g)*s+2.0*c*b;
438  p=s*r;
439  d[i+1]=g+p;
440  g=c*r-b;
441  /* Next loop can be omitted if eigenvectors not wanted */
442  for (int k=1;k<=n;k++)
443  {
444  f=z[k][i+1];
445  z[k][i+1]=s*z[k][i]+c*f;
446  z[k][i]=c*z[k][i]-s*f;
447  }
448  }
449  d[l]=d[l]-p;
450  e[l]=g;
451  e[m]=0.0;
452  }
453  } while (m != l);
454  }
455  return d;
456 }
457 #undef EIGEN_VECTORS
dvar_vector get_eigen_values(const dvar_vector &_ddd, const dvar_vector &_eee)
Eigenvalues and eigenvectors.
Definition: dveigen.cpp:269
dvector eigenvalues(const banded_symmetric_dmatrix &_SS)
Description not yet available.
Definition: dmat28.cpp:411
#define x
Vector of double precision numbers.
Definition: dvector.h:50
int indexmin() const
Get minimum valid index.
Definition: dvector.h:199
df1_two_variable fabs(const df1_two_variable &x)
Definition: df12fun.cpp:891
exitptr ad_exit
Definition: gradstrc.cpp:53
unsigned int colsize() const
Definition: fvar.hpp:2948
void get_eigen(const dvar_vector &d, const dvar_vector &e, const dvar_matrix &z)
Eigenvalues.
Definition: dveigen.cpp:187
dmatrix symmetrize(const dmatrix &matrix)
Author: David Fournier Copyright (c) 2008-2012 Regents of the University of California.
Definition: dmat12.cpp:14
d3_array sqrt(const d3_array &arr3)
Author: David Fournier Copyright (c) 2008-2012 Regents of the University of California.
Definition: d3arr2c.cpp:11
int indexmax() const
Get maximum valid index.
Definition: dvector.h:204
Author: David Fournier Copyright (c) 2008-2012 Regents of the University of California.
Description not yet available.
Definition: fvar.hpp:2819
unsigned int size() const
Get number of elements in array.
Definition: dvector.h:209
void tri_dag(const dvar_matrix &, const dvar_vector &, const dvar_vector &)
Householder transformation for eigenvalue computation.
Definition: dveigen.cpp:54
void rowshift(int min)
Changes the range of valid indices for the rows.
Definition: dmat9.cpp:48
void colshift(int min)
Description not yet available.
Definition: dmat9.cpp:68
unsigned int rowsize() const
Definition: fvar.hpp:2934
#define SIGN(a, b)
Definition: dmat42.cpp:20