For God so loved the world, that He gave His only begotten Son, that all who believe in Him should not perish but have everlasting life
You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

528 lines
18 KiB

  1. <?php
  2. /**
  3. * @package JAMA
  4. *
  5. * For an m-by-n matrix A with m >= n, the singular value decomposition is
  6. * an m-by-n orthogonal matrix U, an n-by-n diagonal matrix S, and
  7. * an n-by-n orthogonal matrix V so that A = U*S*V'.
  8. *
  9. * The singular values, sigma[$k] = S[$k][$k], are ordered so that
  10. * sigma[0] >= sigma[1] >= ... >= sigma[n-1].
  11. *
  12. * The singular value decompostion always exists, so the constructor will
  13. * never fail. The matrix condition number and the effective numerical
  14. * rank can be computed from this decomposition.
  15. *
  16. * @author Paul Meagher
  17. * @license PHP v3.0
  18. * @version 1.1
  19. */
  20. class SingularValueDecomposition
  21. {
  22. /**
  23. * Internal storage of U.
  24. * @var array
  25. */
  26. private $U = array();
  27. /**
  28. * Internal storage of V.
  29. * @var array
  30. */
  31. private $V = array();
  32. /**
  33. * Internal storage of singular values.
  34. * @var array
  35. */
  36. private $s = array();
  37. /**
  38. * Row dimension.
  39. * @var int
  40. */
  41. private $m;
  42. /**
  43. * Column dimension.
  44. * @var int
  45. */
  46. private $n;
  47. /**
  48. * Construct the singular value decomposition
  49. *
  50. * Derived from LINPACK code.
  51. *
  52. * @param $A Rectangular matrix
  53. * @return Structure to access U, S and V.
  54. */
  55. public function __construct($Arg)
  56. {
  57. // Initialize.
  58. $A = $Arg->getArrayCopy();
  59. $this->m = $Arg->getRowDimension();
  60. $this->n = $Arg->getColumnDimension();
  61. $nu = min($this->m, $this->n);
  62. $e = array();
  63. $work = array();
  64. $wantu = true;
  65. $wantv = true;
  66. $nct = min($this->m - 1, $this->n);
  67. $nrt = max(0, min($this->n - 2, $this->m));
  68. // Reduce A to bidiagonal form, storing the diagonal elements
  69. // in s and the super-diagonal elements in e.
  70. for ($k = 0; $k < max($nct, $nrt); ++$k) {
  71. if ($k < $nct) {
  72. // Compute the transformation for the k-th column and
  73. // place the k-th diagonal in s[$k].
  74. // Compute 2-norm of k-th column without under/overflow.
  75. $this->s[$k] = 0;
  76. for ($i = $k; $i < $this->m; ++$i) {
  77. $this->s[$k] = hypo($this->s[$k], $A[$i][$k]);
  78. }
  79. if ($this->s[$k] != 0.0) {
  80. if ($A[$k][$k] < 0.0) {
  81. $this->s[$k] = -$this->s[$k];
  82. }
  83. for ($i = $k; $i < $this->m; ++$i) {
  84. $A[$i][$k] /= $this->s[$k];
  85. }
  86. $A[$k][$k] += 1.0;
  87. }
  88. $this->s[$k] = -$this->s[$k];
  89. }
  90. for ($j = $k + 1; $j < $this->n; ++$j) {
  91. if (($k < $nct) & ($this->s[$k] != 0.0)) {
  92. // Apply the transformation.
  93. $t = 0;
  94. for ($i = $k; $i < $this->m; ++$i) {
  95. $t += $A[$i][$k] * $A[$i][$j];
  96. }
  97. $t = -$t / $A[$k][$k];
  98. for ($i = $k; $i < $this->m; ++$i) {
  99. $A[$i][$j] += $t * $A[$i][$k];
  100. }
  101. // Place the k-th row of A into e for the
  102. // subsequent calculation of the row transformation.
  103. $e[$j] = $A[$k][$j];
  104. }
  105. }
  106. if ($wantu and ($k < $nct)) {
  107. // Place the transformation in U for subsequent back
  108. // multiplication.
  109. for ($i = $k; $i < $this->m; ++$i) {
  110. $this->U[$i][$k] = $A[$i][$k];
  111. }
  112. }
  113. if ($k < $nrt) {
  114. // Compute the k-th row transformation and place the
  115. // k-th super-diagonal in e[$k].
  116. // Compute 2-norm without under/overflow.
  117. $e[$k] = 0;
  118. for ($i = $k + 1; $i < $this->n; ++$i) {
  119. $e[$k] = hypo($e[$k], $e[$i]);
  120. }
  121. if ($e[$k] != 0.0) {
  122. if ($e[$k+1] < 0.0) {
  123. $e[$k] = -$e[$k];
  124. }
  125. for ($i = $k + 1; $i < $this->n; ++$i) {
  126. $e[$i] /= $e[$k];
  127. }
  128. $e[$k+1] += 1.0;
  129. }
  130. $e[$k] = -$e[$k];
  131. if (($k+1 < $this->m) and ($e[$k] != 0.0)) {
  132. // Apply the transformation.
  133. for ($i = $k+1; $i < $this->m; ++$i) {
  134. $work[$i] = 0.0;
  135. }
  136. for ($j = $k+1; $j < $this->n; ++$j) {
  137. for ($i = $k+1; $i < $this->m; ++$i) {
  138. $work[$i] += $e[$j] * $A[$i][$j];
  139. }
  140. }
  141. for ($j = $k + 1; $j < $this->n; ++$j) {
  142. $t = -$e[$j] / $e[$k+1];
  143. for ($i = $k + 1; $i < $this->m; ++$i) {
  144. $A[$i][$j] += $t * $work[$i];
  145. }
  146. }
  147. }
  148. if ($wantv) {
  149. // Place the transformation in V for subsequent
  150. // back multiplication.
  151. for ($i = $k + 1; $i < $this->n; ++$i) {
  152. $this->V[$i][$k] = $e[$i];
  153. }
  154. }
  155. }
  156. }
  157. // Set up the final bidiagonal matrix or order p.
  158. $p = min($this->n, $this->m + 1);
  159. if ($nct < $this->n) {
  160. $this->s[$nct] = $A[$nct][$nct];
  161. }
  162. if ($this->m < $p) {
  163. $this->s[$p-1] = 0.0;
  164. }
  165. if ($nrt + 1 < $p) {
  166. $e[$nrt] = $A[$nrt][$p-1];
  167. }
  168. $e[$p-1] = 0.0;
  169. // If required, generate U.
  170. if ($wantu) {
  171. for ($j = $nct; $j < $nu; ++$j) {
  172. for ($i = 0; $i < $this->m; ++$i) {
  173. $this->U[$i][$j] = 0.0;
  174. }
  175. $this->U[$j][$j] = 1.0;
  176. }
  177. for ($k = $nct - 1; $k >= 0; --$k) {
  178. if ($this->s[$k] != 0.0) {
  179. for ($j = $k + 1; $j < $nu; ++$j) {
  180. $t = 0;
  181. for ($i = $k; $i < $this->m; ++$i) {
  182. $t += $this->U[$i][$k] * $this->U[$i][$j];
  183. }
  184. $t = -$t / $this->U[$k][$k];
  185. for ($i = $k; $i < $this->m; ++$i) {
  186. $this->U[$i][$j] += $t * $this->U[$i][$k];
  187. }
  188. }
  189. for ($i = $k; $i < $this->m; ++$i) {
  190. $this->U[$i][$k] = -$this->U[$i][$k];
  191. }
  192. $this->U[$k][$k] = 1.0 + $this->U[$k][$k];
  193. for ($i = 0; $i < $k - 1; ++$i) {
  194. $this->U[$i][$k] = 0.0;
  195. }
  196. } else {
  197. for ($i = 0; $i < $this->m; ++$i) {
  198. $this->U[$i][$k] = 0.0;
  199. }
  200. $this->U[$k][$k] = 1.0;
  201. }
  202. }
  203. }
  204. // If required, generate V.
  205. if ($wantv) {
  206. for ($k = $this->n - 1; $k >= 0; --$k) {
  207. if (($k < $nrt) and ($e[$k] != 0.0)) {
  208. for ($j = $k + 1; $j < $nu; ++$j) {
  209. $t = 0;
  210. for ($i = $k + 1; $i < $this->n; ++$i) {
  211. $t += $this->V[$i][$k]* $this->V[$i][$j];
  212. }
  213. $t = -$t / $this->V[$k+1][$k];
  214. for ($i = $k + 1; $i < $this->n; ++$i) {
  215. $this->V[$i][$j] += $t * $this->V[$i][$k];
  216. }
  217. }
  218. }
  219. for ($i = 0; $i < $this->n; ++$i) {
  220. $this->V[$i][$k] = 0.0;
  221. }
  222. $this->V[$k][$k] = 1.0;
  223. }
  224. }
  225. // Main iteration loop for the singular values.
  226. $pp = $p - 1;
  227. $iter = 0;
  228. $eps = pow(2.0, -52.0);
  229. while ($p > 0) {
  230. // Here is where a test for too many iterations would go.
  231. // This section of the program inspects for negligible
  232. // elements in the s and e arrays. On completion the
  233. // variables kase and k are set as follows:
  234. // kase = 1 if s(p) and e[k-1] are negligible and k<p
  235. // kase = 2 if s(k) is negligible and k<p
  236. // kase = 3 if e[k-1] is negligible, k<p, and
  237. // s(k), ..., s(p) are not negligible (qr step).
  238. // kase = 4 if e(p-1) is negligible (convergence).
  239. for ($k = $p - 2; $k >= -1; --$k) {
  240. if ($k == -1) {
  241. break;
  242. }
  243. if (abs($e[$k]) <= $eps * (abs($this->s[$k]) + abs($this->s[$k+1]))) {
  244. $e[$k] = 0.0;
  245. break;
  246. }
  247. }
  248. if ($k == $p - 2) {
  249. $kase = 4;
  250. } else {
  251. for ($ks = $p - 1; $ks >= $k; --$ks) {
  252. if ($ks == $k) {
  253. break;
  254. }
  255. $t = ($ks != $p ? abs($e[$ks]) : 0.) + ($ks != $k + 1 ? abs($e[$ks-1]) : 0.);
  256. if (abs($this->s[$ks]) <= $eps * $t) {
  257. $this->s[$ks] = 0.0;
  258. break;
  259. }
  260. }
  261. if ($ks == $k) {
  262. $kase = 3;
  263. } elseif ($ks == $p-1) {
  264. $kase = 1;
  265. } else {
  266. $kase = 2;
  267. $k = $ks;
  268. }
  269. }
  270. ++$k;
  271. // Perform the task indicated by kase.
  272. switch ($kase) {
  273. // Deflate negligible s(p).
  274. case 1:
  275. $f = $e[$p-2];
  276. $e[$p-2] = 0.0;
  277. for ($j = $p - 2; $j >= $k; --$j) {
  278. $t = hypo($this->s[$j], $f);
  279. $cs = $this->s[$j] / $t;
  280. $sn = $f / $t;
  281. $this->s[$j] = $t;
  282. if ($j != $k) {
  283. $f = -$sn * $e[$j-1];
  284. $e[$j-1] = $cs * $e[$j-1];
  285. }
  286. if ($wantv) {
  287. for ($i = 0; $i < $this->n; ++$i) {
  288. $t = $cs * $this->V[$i][$j] + $sn * $this->V[$i][$p-1];
  289. $this->V[$i][$p-1] = -$sn * $this->V[$i][$j] + $cs * $this->V[$i][$p-1];
  290. $this->V[$i][$j] = $t;
  291. }
  292. }
  293. }
  294. break;
  295. // Split at negligible s(k).
  296. case 2:
  297. $f = $e[$k-1];
  298. $e[$k-1] = 0.0;
  299. for ($j = $k; $j < $p; ++$j) {
  300. $t = hypo($this->s[$j], $f);
  301. $cs = $this->s[$j] / $t;
  302. $sn = $f / $t;
  303. $this->s[$j] = $t;
  304. $f = -$sn * $e[$j];
  305. $e[$j] = $cs * $e[$j];
  306. if ($wantu) {
  307. for ($i = 0; $i < $this->m; ++$i) {
  308. $t = $cs * $this->U[$i][$j] + $sn * $this->U[$i][$k-1];
  309. $this->U[$i][$k-1] = -$sn * $this->U[$i][$j] + $cs * $this->U[$i][$k-1];
  310. $this->U[$i][$j] = $t;
  311. }
  312. }
  313. }
  314. break;
  315. // Perform one qr step.
  316. case 3:
  317. // Calculate the shift.
  318. $scale = max(max(max(max(abs($this->s[$p-1]), abs($this->s[$p-2])), abs($e[$p-2])), abs($this->s[$k])), abs($e[$k]));
  319. $sp = $this->s[$p-1] / $scale;
  320. $spm1 = $this->s[$p-2] / $scale;
  321. $epm1 = $e[$p-2] / $scale;
  322. $sk = $this->s[$k] / $scale;
  323. $ek = $e[$k] / $scale;
  324. $b = (($spm1 + $sp) * ($spm1 - $sp) + $epm1 * $epm1) / 2.0;
  325. $c = ($sp * $epm1) * ($sp * $epm1);
  326. $shift = 0.0;
  327. if (($b != 0.0) || ($c != 0.0)) {
  328. $shift = sqrt($b * $b + $c);
  329. if ($b < 0.0) {
  330. $shift = -$shift;
  331. }
  332. $shift = $c / ($b + $shift);
  333. }
  334. $f = ($sk + $sp) * ($sk - $sp) + $shift;
  335. $g = $sk * $ek;
  336. // Chase zeros.
  337. for ($j = $k; $j < $p-1; ++$j) {
  338. $t = hypo($f, $g);
  339. $cs = $f/$t;
  340. $sn = $g/$t;
  341. if ($j != $k) {
  342. $e[$j-1] = $t;
  343. }
  344. $f = $cs * $this->s[$j] + $sn * $e[$j];
  345. $e[$j] = $cs * $e[$j] - $sn * $this->s[$j];
  346. $g = $sn * $this->s[$j+1];
  347. $this->s[$j+1] = $cs * $this->s[$j+1];
  348. if ($wantv) {
  349. for ($i = 0; $i < $this->n; ++$i) {
  350. $t = $cs * $this->V[$i][$j] + $sn * $this->V[$i][$j+1];
  351. $this->V[$i][$j+1] = -$sn * $this->V[$i][$j] + $cs * $this->V[$i][$j+1];
  352. $this->V[$i][$j] = $t;
  353. }
  354. }
  355. $t = hypo($f, $g);
  356. $cs = $f/$t;
  357. $sn = $g/$t;
  358. $this->s[$j] = $t;
  359. $f = $cs * $e[$j] + $sn * $this->s[$j+1];
  360. $this->s[$j+1] = -$sn * $e[$j] + $cs * $this->s[$j+1];
  361. $g = $sn * $e[$j+1];
  362. $e[$j+1] = $cs * $e[$j+1];
  363. if ($wantu && ($j < $this->m - 1)) {
  364. for ($i = 0; $i < $this->m; ++$i) {
  365. $t = $cs * $this->U[$i][$j] + $sn * $this->U[$i][$j+1];
  366. $this->U[$i][$j+1] = -$sn * $this->U[$i][$j] + $cs * $this->U[$i][$j+1];
  367. $this->U[$i][$j] = $t;
  368. }
  369. }
  370. }
  371. $e[$p-2] = $f;
  372. $iter = $iter + 1;
  373. break;
  374. // Convergence.
  375. case 4:
  376. // Make the singular values positive.
  377. if ($this->s[$k] <= 0.0) {
  378. $this->s[$k] = ($this->s[$k] < 0.0 ? -$this->s[$k] : 0.0);
  379. if ($wantv) {
  380. for ($i = 0; $i <= $pp; ++$i) {
  381. $this->V[$i][$k] = -$this->V[$i][$k];
  382. }
  383. }
  384. }
  385. // Order the singular values.
  386. while ($k < $pp) {
  387. if ($this->s[$k] >= $this->s[$k+1]) {
  388. break;
  389. }
  390. $t = $this->s[$k];
  391. $this->s[$k] = $this->s[$k+1];
  392. $this->s[$k+1] = $t;
  393. if ($wantv and ($k < $this->n - 1)) {
  394. for ($i = 0; $i < $this->n; ++$i) {
  395. $t = $this->V[$i][$k+1];
  396. $this->V[$i][$k+1] = $this->V[$i][$k];
  397. $this->V[$i][$k] = $t;
  398. }
  399. }
  400. if ($wantu and ($k < $this->m-1)) {
  401. for ($i = 0; $i < $this->m; ++$i) {
  402. $t = $this->U[$i][$k+1];
  403. $this->U[$i][$k+1] = $this->U[$i][$k];
  404. $this->U[$i][$k] = $t;
  405. }
  406. }
  407. ++$k;
  408. }
  409. $iter = 0;
  410. --$p;
  411. break;
  412. } // end switch
  413. } // end while
  414. } // end constructor
  415. /**
  416. * Return the left singular vectors
  417. *
  418. * @access public
  419. * @return U
  420. */
  421. public function getU()
  422. {
  423. return new Matrix($this->U, $this->m, min($this->m + 1, $this->n));
  424. }
  425. /**
  426. * Return the right singular vectors
  427. *
  428. * @access public
  429. * @return V
  430. */
  431. public function getV()
  432. {
  433. return new Matrix($this->V);
  434. }
  435. /**
  436. * Return the one-dimensional array of singular values
  437. *
  438. * @access public
  439. * @return diagonal of S.
  440. */
  441. public function getSingularValues()
  442. {
  443. return $this->s;
  444. }
  445. /**
  446. * Return the diagonal matrix of singular values
  447. *
  448. * @access public
  449. * @return S
  450. */
  451. public function getS()
  452. {
  453. for ($i = 0; $i < $this->n; ++$i) {
  454. for ($j = 0; $j < $this->n; ++$j) {
  455. $S[$i][$j] = 0.0;
  456. }
  457. $S[$i][$i] = $this->s[$i];
  458. }
  459. return new Matrix($S);
  460. }
  461. /**
  462. * Two norm
  463. *
  464. * @access public
  465. * @return max(S)
  466. */
  467. public function norm2()
  468. {
  469. return $this->s[0];
  470. }
  471. /**
  472. * Two norm condition number
  473. *
  474. * @access public
  475. * @return max(S)/min(S)
  476. */
  477. public function cond()
  478. {
  479. return $this->s[0] / $this->s[min($this->m, $this->n) - 1];
  480. }
  481. /**
  482. * Effective numerical matrix rank
  483. *
  484. * @access public
  485. * @return Number of nonnegligible singular values.
  486. */
  487. public function rank()
  488. {
  489. $eps = pow(2.0, -52.0);
  490. $tol = max($this->m, $this->n) * $this->s[0] * $eps;
  491. $r = 0;
  492. for ($i = 0; $i < count($this->s); ++$i) {
  493. if ($this->s[$i] > $tol) {
  494. ++$r;
  495. }
  496. }
  497. return $r;
  498. }
  499. }