1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23 package org.hipparchus.linear;
24
25 import org.hipparchus.exception.LocalizedCoreFormats;
26 import org.hipparchus.exception.MathIllegalArgumentException;
27 import org.hipparchus.exception.MathIllegalStateException;
28 import org.hipparchus.exception.MathRuntimeException;
29 import org.hipparchus.util.FastMath;
30 import org.hipparchus.util.Precision;
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 public class EigenDecompositionSymmetric {
73
74
75 public static final double DEFAULT_EPSILON = 1e-12;
76
77
78 private static final byte MAX_ITER = 30;
79
80
81 private final double epsilon;
82
83
84 private double[] eigenvalues;
85
86
87 private ArrayRealVector[] eigenvectors;
88
89
90 private RealMatrix cachedV;
91
92
93 private DiagonalMatrix cachedD;
94
95
96 private RealMatrix cachedVt;
97
98
99
100
101
102
103
104
105
106
107
108
109 public EigenDecompositionSymmetric(final RealMatrix matrix) {
110 this(matrix, DEFAULT_EPSILON, true);
111 }
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126 public EigenDecompositionSymmetric(final RealMatrix matrix,
127 final double epsilon, final boolean decreasing)
128 throws MathRuntimeException {
129
130 this.epsilon = epsilon;
131 MatrixUtils.checkSymmetric(matrix, epsilon);
132
133
134 final TriDiagonalTransformer transformer = new TriDiagonalTransformer(matrix);
135
136 findEigenVectors(transformer.getMainDiagonalRef(),
137 transformer.getSecondaryDiagonalRef(),
138 transformer.getQ().getData(),
139 decreasing);
140
141 }
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156 public EigenDecompositionSymmetric(final double[] main, final double[] secondary) {
157 this(main, secondary, DEFAULT_EPSILON, true);
158 }
159
160
161
162
163
164
165
166
167
168
169
170
171 public EigenDecompositionSymmetric(final double[] main, final double[] secondary,
172 final double epsilon, final boolean decreasing) {
173 this.epsilon = epsilon;
174 final int size = main.length;
175 final double[][] z = new double[size][size];
176 for (int i = 0; i < size; i++) {
177 z[i][i] = 1.0;
178 }
179 findEigenVectors(main.clone(), secondary.clone(), z, decreasing);
180 }
181
182
183
184
185
186
187
188
189
190
191
192 public RealMatrix getV() {
193
194 if (cachedV == null) {
195 final int m = eigenvectors.length;
196 cachedV = MatrixUtils.createRealMatrix(m, m);
197 for (int k = 0; k < m; ++k) {
198 cachedV.setColumnVector(k, eigenvectors[k]);
199 }
200 }
201
202 return cachedV;
203 }
204
205
206
207
208
209
210
211
212 public DiagonalMatrix getD() {
213
214 if (cachedD == null) {
215
216 cachedD = new DiagonalMatrix(eigenvalues);
217 }
218
219 return cachedD;
220
221 }
222
223
224
225
226
227
228 public double getEpsilon() { return epsilon; }
229
230
231
232
233
234
235
236
237
238
239
240 public RealMatrix getVT() {
241
242 if (cachedVt == null) {
243 final int m = eigenvectors.length;
244 cachedVt = MatrixUtils.createRealMatrix(m, m);
245 for (int k = 0; k < m; ++k) {
246 cachedVt.setRowVector(k, eigenvectors[k]);
247 }
248 }
249
250
251 return cachedVt;
252 }
253
254
255
256
257
258
259
260
261
262 public double[] getEigenvalues() {
263 return eigenvalues.clone();
264 }
265
266
267
268
269
270
271
272
273
274
275
276 public double getEigenvalue(final int i) {
277 return eigenvalues[i];
278 }
279
280
281
282
283
284
285
286
287
288
289
290 public RealVector getEigenvector(final int i) {
291 return eigenvectors[i].copy();
292 }
293
294
295
296
297
298
299 public double getDeterminant() {
300 double determinant = 1;
301 for (int i = 0; i < eigenvalues.length; ++i) {
302 determinant *= eigenvalues[i];
303 }
304 return determinant;
305 }
306
307
308
309
310
311
312
313
314
315 public RealMatrix getSquareRoot() {
316
317 final double[] sqrtEigenValues = new double[eigenvalues.length];
318 for (int i = 0; i < eigenvalues.length; i++) {
319 final double eigen = eigenvalues[i];
320 if (eigen <= 0) {
321 throw new MathRuntimeException(LocalizedCoreFormats.UNSUPPORTED_OPERATION);
322 }
323 sqrtEigenValues[i] = FastMath.sqrt(eigen);
324 }
325 final RealMatrix sqrtEigen = MatrixUtils.createRealDiagonalMatrix(sqrtEigenValues);
326 final RealMatrix v = getV();
327 final RealMatrix vT = getVT();
328
329 return v.multiply(sqrtEigen).multiply(vT);
330
331 }
332
333
334
335
336 public DecompositionSolver getSolver() {
337 return new Solver();
338 }
339
340
341 private class Solver implements DecompositionSolver {
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356 @Override
357 public RealVector solve(final RealVector b) {
358 if (!isNonSingular()) {
359 throw new MathIllegalArgumentException(LocalizedCoreFormats.SINGULAR_MATRIX);
360 }
361
362 final int m = eigenvalues.length;
363 if (b.getDimension() != m) {
364 throw new MathIllegalArgumentException(LocalizedCoreFormats.DIMENSIONS_MISMATCH,
365 b.getDimension(), m);
366 }
367
368 final double[] bp = new double[m];
369 for (int i = 0; i < m; ++i) {
370 final ArrayRealVector v = eigenvectors[i];
371 final double[] vData = v.getDataRef();
372 final double s = v.dotProduct(b) / eigenvalues[i];
373 for (int j = 0; j < m; ++j) {
374 bp[j] += s * vData[j];
375 }
376 }
377
378 return new ArrayRealVector(bp, false);
379 }
380
381
382 @Override
383 public RealMatrix solve(RealMatrix b) {
384
385 if (!isNonSingular()) {
386 throw new MathIllegalArgumentException(LocalizedCoreFormats.SINGULAR_MATRIX);
387 }
388
389 final int m = eigenvalues.length;
390 if (b.getRowDimension() != m) {
391 throw new MathIllegalArgumentException(LocalizedCoreFormats.DIMENSIONS_MISMATCH,
392 b.getRowDimension(), m);
393 }
394
395 final int nColB = b.getColumnDimension();
396 final double[][] bp = new double[m][nColB];
397 final double[] tmpCol = new double[m];
398 for (int k = 0; k < nColB; ++k) {
399 for (int i = 0; i < m; ++i) {
400 tmpCol[i] = b.getEntry(i, k);
401 bp[i][k] = 0;
402 }
403 for (int i = 0; i < m; ++i) {
404 final ArrayRealVector v = eigenvectors[i];
405 final double[] vData = v.getDataRef();
406 double s = 0;
407 for (int j = 0; j < m; ++j) {
408 s += v.getEntry(j) * tmpCol[j];
409 }
410 s /= eigenvalues[i];
411 for (int j = 0; j < m; ++j) {
412 bp[j][k] += s * vData[j];
413 }
414 }
415 }
416
417 return new Array2DRowRealMatrix(bp, false);
418
419 }
420
421
422
423
424
425
426 @Override
427 public boolean isNonSingular() {
428 double largestEigenvalueNorm = 0.0;
429
430
431 for (int i = 0; i < eigenvalues.length; ++i) {
432 largestEigenvalueNorm = FastMath.max(largestEigenvalueNorm, FastMath.abs(eigenvalues[i]));
433 }
434
435 if (largestEigenvalueNorm == 0.0) {
436 return false;
437 }
438 for (int i = 0; i < eigenvalues.length; ++i) {
439
440
441 if (Precision.equals(FastMath.abs(eigenvalues[i]) / largestEigenvalueNorm, 0, epsilon)) {
442 return false;
443 }
444 }
445 return true;
446 }
447
448
449
450
451
452
453
454 @Override
455 public RealMatrix getInverse() {
456 if (!isNonSingular()) {
457 throw new MathIllegalArgumentException(LocalizedCoreFormats.SINGULAR_MATRIX);
458 }
459
460 final int m = eigenvalues.length;
461 final double[][] invData = new double[m][m];
462
463 for (int i = 0; i < m; ++i) {
464 final double[] invI = invData[i];
465 for (int j = 0; j < m; ++j) {
466 double invIJ = 0;
467 for (int k = 0; k < m; ++k) {
468 final double[] vK = eigenvectors[k].getDataRef();
469 invIJ += vK[i] * vK[j] / eigenvalues[k];
470 }
471 invI[j] = invIJ;
472 }
473 }
474 return MatrixUtils.createRealMatrix(invData);
475 }
476
477
478 @Override
479 public int getRowDimension() {
480 return eigenvalues.length;
481 }
482
483
484 @Override
485 public int getColumnDimension() {
486 return eigenvalues.length;
487 }
488
489 }
490
491
492
493
494
495
496
497
498
499 private void findEigenVectors(final double[] main, final double[] secondary,
500 final double[][] householderMatrix, final boolean decreasing) {
501 final double[][]z = householderMatrix.clone();
502 final int n = main.length;
503 eigenvalues = new double[n];
504 final double[] e = new double[n];
505 for (int i = 0; i < n - 1; i++) {
506 eigenvalues[i] = main[i];
507 e[i] = secondary[i];
508 }
509 eigenvalues[n - 1] = main[n - 1];
510 e[n - 1] = 0;
511
512
513 double maxAbsoluteValue = 0;
514 for (int i = 0; i < n; i++) {
515 if (FastMath.abs(eigenvalues[i]) > maxAbsoluteValue) {
516 maxAbsoluteValue = FastMath.abs(eigenvalues[i]);
517 }
518 if (FastMath.abs(e[i]) > maxAbsoluteValue) {
519 maxAbsoluteValue = FastMath.abs(e[i]);
520 }
521 }
522
523 if (maxAbsoluteValue != 0) {
524 for (int i=0; i < n; i++) {
525 if (FastMath.abs(eigenvalues[i]) <= Precision.EPSILON * maxAbsoluteValue) {
526 eigenvalues[i] = 0;
527 }
528 if (FastMath.abs(e[i]) <= Precision.EPSILON * maxAbsoluteValue) {
529 e[i]=0;
530 }
531 }
532 }
533
534 for (int j = 0; j < n; j++) {
535 int its = 0;
536 int m;
537 do {
538 for (m = j; m < n - 1; m++) {
539 double delta = FastMath.abs(eigenvalues[m]) +
540 FastMath.abs(eigenvalues[m + 1]);
541 if (FastMath.abs(e[m]) + delta == delta) {
542 break;
543 }
544 }
545 if (m != j) {
546 if (its == MAX_ITER) {
547 throw new MathIllegalStateException(LocalizedCoreFormats.CONVERGENCE_FAILED,
548 MAX_ITER);
549 }
550 its++;
551 double q = (eigenvalues[j + 1] - eigenvalues[j]) / (2 * e[j]);
552 double t = FastMath.sqrt(1 + q * q);
553 if (q < 0.0) {
554 q = eigenvalues[m] - eigenvalues[j] + e[j] / (q - t);
555 } else {
556 q = eigenvalues[m] - eigenvalues[j] + e[j] / (q + t);
557 }
558 double u = 0.0;
559 double s = 1.0;
560 double c = 1.0;
561 int i;
562 for (i = m - 1; i >= j; i--) {
563 double p = s * e[i];
564 double h = c * e[i];
565 if (FastMath.abs(p) >= FastMath.abs(q)) {
566 c = q / p;
567 t = FastMath.sqrt(c * c + 1.0);
568 e[i + 1] = p * t;
569 s = 1.0 / t;
570 c *= s;
571 } else {
572 s = p / q;
573 t = FastMath.sqrt(s * s + 1.0);
574 e[i + 1] = q * t;
575 c = 1.0 / t;
576 s *= c;
577 }
578 if (e[i + 1] == 0.0) {
579 eigenvalues[i + 1] -= u;
580 e[m] = 0.0;
581 break;
582 }
583 q = eigenvalues[i + 1] - u;
584 t = (eigenvalues[i] - q) * s + 2.0 * c * h;
585 u = s * t;
586 eigenvalues[i + 1] = q + u;
587 q = c * t - h;
588 for (int ia = 0; ia < n; ia++) {
589 p = z[ia][i + 1];
590 z[ia][i + 1] = s * z[ia][i] + c * p;
591 z[ia][i] = c * z[ia][i] - s * p;
592 }
593 }
594 if (t == 0.0 && i >= j) {
595 continue;
596 }
597 eigenvalues[j] -= u;
598 e[j] = q;
599 e[m] = 0.0;
600 }
601 } while (m != j);
602 }
603
604
605 for (int i = 0; i < n; i++) {
606 int k = i;
607 double p = eigenvalues[i];
608 for (int j = i + 1; j < n; j++) {
609 if (eigenvalues[j] > p == decreasing) {
610 k = j;
611 p = eigenvalues[j];
612 }
613 }
614 if (k != i) {
615 eigenvalues[k] = eigenvalues[i];
616 eigenvalues[i] = p;
617 for (int j = 0; j < n; j++) {
618 p = z[j][i];
619 z[j][i] = z[j][k];
620 z[j][k] = p;
621 }
622 }
623 }
624
625
626 maxAbsoluteValue = 0;
627 for (int i = 0; i < n; i++) {
628 if (FastMath.abs(eigenvalues[i]) > maxAbsoluteValue) {
629 maxAbsoluteValue = FastMath.abs(eigenvalues[i]);
630 }
631 }
632
633 if (maxAbsoluteValue != 0.0) {
634 for (int i=0; i < n; i++) {
635 if (FastMath.abs(eigenvalues[i]) < Precision.EPSILON * maxAbsoluteValue) {
636 eigenvalues[i] = 0;
637 }
638 }
639 }
640 eigenvectors = new ArrayRealVector[n];
641 for (int i = 0; i < n; i++) {
642 eigenvectors[i] = new ArrayRealVector(n);
643 for (int j = 0; j < n; j++) {
644 eigenvectors[i].setEntry(j, z[j][i]);
645 }
646 }
647 }
648
649 }