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
|
// File: math_EigenValuesSearcher.cxx
// Created: Thu Dec 15 17:58:17 2005
// Author: Julia GERASIMOVA
// <jgv@clubox>
#include <math_EigenValuesSearcher.ixx>
//==========================================================================
//function : pythag
// Computation of sqrt(x*x + y*y).
//==========================================================================
static inline Standard_Real pythag(const Standard_Real x,
const Standard_Real y)
{
return Sqrt(x*x + y*y);
}
math_EigenValuesSearcher::math_EigenValuesSearcher(const TColStd_Array1OfReal& Diagonal,
const TColStd_Array1OfReal& Subdiagonal)
{
myIsDone = Standard_False;
Standard_Integer n = Diagonal.Length();
if (Subdiagonal.Length() != n)
Standard_Failure::Raise("math_EigenValuesSearcher : dimension mismatch");
myDiagonal = new TColStd_HArray1OfReal(1, n);
myDiagonal->ChangeArray1() = Diagonal;
mySubdiagonal = new TColStd_HArray1OfReal(1, n);
mySubdiagonal->ChangeArray1() = Subdiagonal;
myN = n;
myEigenValues = new TColStd_HArray1OfReal(1, n);
myEigenVectors = new TColStd_HArray2OfReal(1, n, 1, n);
Standard_Real* d = new Standard_Real [n+1];
Standard_Real* e = new Standard_Real [n+1];
Standard_Real** z = new Standard_Real* [n+1];
Standard_Integer i, j;
for (i = 1; i <= n; i++)
z[i] = new Standard_Real [n+1];
for (i = 1; i <= n; i++)
d[i] = myDiagonal->Value(i);
for (i = 2; i <= n; i++)
e[i] = mySubdiagonal->Value(i);
for (i = 1; i <= n; i++)
for (j = 1; j <= n; j++)
z[i][j] = (i == j)? 1. : 0.;
Standard_Boolean result;
Standard_Integer m;
Standard_Integer l;
Standard_Integer iter;
//Standard_Integer i;
Standard_Integer k;
Standard_Real s;
Standard_Real r;
Standard_Real p;
Standard_Real g;
Standard_Real f;
Standard_Real dd;
Standard_Real c;
Standard_Real b;
result = Standard_True;
if (n != 1)
{
// Shift e.
for (i = 2; i <= n; i++)
e[i - 1] = e[i];
e[n] = 0.0;
for (l = 1; l <= n; l++) {
iter = 0;
do {
for (m = l; m <= n-1; m++) {
dd = Abs(d[m]) + Abs(d[m + 1]);
if (Abs(e[m]) + dd == dd)
break;
}
if (m != l) {
if (iter++ == 30) {
result = Standard_False;
break; //return result;
}
g = (d[l + 1] - d[l])/(2.*e[l]);
r = pythag(1., g);
if (g < 0)
g = d[m] - d[l] + e[l]/(g - r);
else
g = d[m] - d[l] + e[l]/(g + r);
s = 1.;
c = 1.;
p = 0.;
for (i = m - 1; i >= l; i--) {
f = s*e[i];
b = c*e[i];
r = pythag (f, g);
e[i + 1] = r;
if (r == 0.) {
d[i + 1] -= p;
e[m] = 0.;
break;
}
s = f/r;
c = g/r;
g = d[i + 1] - p;
r = (d[i] - g)*s + 2.0*c*b;
p = s*r;
d[i + 1] = g + p;
g = c*r - b;
for (k = 1; k <= n; k++) {
f = z[k][i + 1];
z[k][i + 1] = s*z[k][i] + c*f;
z[k][i] = c*z[k][i] - s*f;
}
}
if (r == 0 && i >= 1)
continue;
d[l] -= p;
e[l] = g;
e[m] = 0.;
}
}
while (m != l);
if (result == Standard_False)
break;
} //end of for (l = 1; l <= n; l++)
} //end of if (n != 1)
if (result)
{
for (i = 1; i <= n; i++)
myEigenValues->ChangeValue(i) = d[i];
for (i = 1; i <= n; i++)
for (j = 1; j <= n; j++)
myEigenVectors->ChangeValue(i, j) = z[i][j];
}
myIsDone = result;
delete [] d;
delete [] e;
for (i = 1; i <= n; i++)
delete [] z[i];
delete [] z;
}
Standard_Boolean math_EigenValuesSearcher::IsDone() const
{
return myIsDone;
}
Standard_Integer math_EigenValuesSearcher::Dimension() const
{
return myN;
}
Standard_Real math_EigenValuesSearcher::EigenValue(const Standard_Integer Index) const
{
return myEigenValues->Value(Index);
}
math_Vector math_EigenValuesSearcher::EigenVector(const Standard_Integer Index) const
{
math_Vector theVector(1, myN);
Standard_Integer i;
for (i = 1; i <= myN; i++)
theVector(i) = myEigenVectors->Value(i, Index);
return theVector;
}
|