pca.m 8.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312
  1. function [U,S,V] = pca(A,k,its,l)
  2. %PCA Low-rank approximation in SVD form.
  3. %
  4. %
  5. % [U,S,V] = PCA(A) constructs a nearly optimal rank-6 approximation
  6. % USV' to A, using 2 full iterations of a block Lanczos method
  7. % of block size 6+2=8, started with an n x 8 random matrix,
  8. % when A is m x n; the ref. below explains "nearly optimal."
  9. % The smallest dimension of A must be >= 6 when A is
  10. % the only input to PCA.
  11. %
  12. % [U,S,V] = PCA(A,k) constructs a nearly optimal rank-k approximation
  13. % USV' to A, using 2 full iterations of a block Lanczos method
  14. % of block size k+2, started with an n x (k+2) random matrix,
  15. % when A is m x n; the ref. below explains "nearly optimal."
  16. % k must be a positive integer <= the smallest dimension of A.
  17. %
  18. % [U,S,V] = PCA(A,k,its) constructs a nearly optimal rank-k approx. USV'
  19. % to A, using its full iterations of a block Lanczos method
  20. % of block size k+2, started with an n x (k+2) random matrix,
  21. % when A is m x n; the ref. below explains "nearly optimal."
  22. % k must be a positive integer <= the smallest dimension of A,
  23. % and its must be a nonnegative integer.
  24. %
  25. % [U,S,V] = PCA(A,k,its,l) constructs a nearly optimal rank-k approx.
  26. % USV' to A, using its full iterates of a block Lanczos method
  27. % of block size l, started with an n x l random matrix,
  28. % when A is m x n; the ref. below explains "nearly optimal."
  29. % k must be a positive integer <= the smallest dimension of A,
  30. % its must be a nonnegative integer,
  31. % and l must be a positive integer >= k.
  32. %
  33. %
  34. % The low-rank approximation USV' is in the form of an SVD in the sense
  35. % that the columns of U are orthonormal, as are the columns of V,
  36. % the entries of S are all nonnegative, and the only nonzero entries
  37. % of S appear in non-increasing order on its diagonal.
  38. % U is m x k, V is n x k, and S is k x k, when A is m x n.
  39. %
  40. % Increasing its or l improves the accuracy of the approximation USV'
  41. % to A; the ref. below describes how the accuracy depends on its and l.
  42. %
  43. %
  44. % Note: PCA invokes RAND. To obtain repeatable results,
  45. % invoke RAND('seed',j) with a fixed integer j before invoking PCA.
  46. %
  47. % Note: PCA currently requires the user to center and normalize the rows
  48. % or columns of the input matrix A before invoking PCA (if such
  49. % is desired).
  50. %
  51. % Note: Hoping to improve its accuracy, PCA takes the adjoint
  52. % of the input matrix A if A has more columns than rows. However,
  53. % the code will produce a correct result without taking an adjoint.
  54. % If taking the adjoint turns out to be costly in a particular
  55. % application, then the user might consider deleting the lines
  56. % below that take the adjoint (assuming its = 0).
  57. %
  58. % Note: The user may ascertain the accuracy of the approximation USV'
  59. % to A by invoking DIFFSNORM(A,U,S,V).
  60. %
  61. %
  62. % inputs (the first is required):
  63. % A -- matrix being approximated
  64. % k -- rank of the approximation being constructed;
  65. % k must be a positive integer <= the smallest dimension of A,
  66. % and defaults to 6
  67. % its -- number of full iterations of a block Lanczos method to conduct;
  68. % its must be a nonnegative integer, and defaults to 2
  69. % l -- block size of the block Lanczos iterations;
  70. % l must be a positive integer >= k, and defaults to k+2
  71. %
  72. % outputs (all three are required):
  73. % U -- m x k matrix in the rank-k approximation USV' to A,
  74. % where A is m x n; the columns of U are orthonormal
  75. % S -- k x k matrix in the rank-k approximation USV' to A,
  76. % where A is m x n; the entries of S are all nonnegative,
  77. % and its only nonzero entries appear in nonincreasing order
  78. % on the diagonal
  79. % V -- n x k matrix in the rank-k approximation USV' to A,
  80. % where A is m x n; the columns of V are orthonormal
  81. %
  82. %
  83. % Example:
  84. % A = rand(1000,2)*rand(2,1000);
  85. % A = A/normest(A);
  86. % [U,S,V] = pca(A,2,0);
  87. % diffsnorm(A,U,S,V)
  88. %
  89. % This code snippet produces a rank-2 approximation USV' to A such that
  90. % the columns of U are orthonormal, as are the columns of V, and
  91. % the entries of S are all nonnegative and are zero off the diagonal.
  92. % diffsnorm(A,U,S,V) outputs an estimate of the spectral norm
  93. % of A-USV', which should be close to the machine precision.
  94. %
  95. %
  96. % Reference:
  97. % Vladimir Rokhlin, Arthur Szlam, and Mark Tygert,
  98. % A randomized algorithm for principal component analysis,
  99. % arXiv:0809.2274v1 [stat.CO], 2008 (available at http://arxiv.org).
  100. %
  101. %
  102. % See also PCACOV, PRINCOMP, SVDS.
  103. %
  104. % Copyright 2008 Mark Tygert.
  105. %
  106. % Check the number of inputs.
  107. %
  108. if(nargin < 1)
  109. error('MATLAB:pca:TooFewIn',...
  110. 'There must be at least 1 input.')
  111. end
  112. if(nargin > 4)
  113. error('MATLAB:pca:TooManyIn',...
  114. 'There must be at most 4 inputs.')
  115. end
  116. %
  117. % Check the number of outputs.
  118. %
  119. if(nargout ~= 3)
  120. error('MATLAB:pca:WrongNumOut',...
  121. 'There must be exactly 3 outputs.')
  122. end
  123. %
  124. % Set the inputs k, its, and l to default values, if necessary.
  125. %
  126. if(nargin == 1)
  127. k = 6;
  128. its = 2;
  129. l = k+2;
  130. end
  131. if(nargin == 2)
  132. its = 2;
  133. l = k+2;
  134. end
  135. if(nargin == 3)
  136. l = k+2;
  137. end
  138. %
  139. % Check the first input argument.
  140. %
  141. if(~isfloat(A))
  142. error('MATLAB:pca:In1NotFloat',...
  143. 'Input 1 must be a floating-point matrix.')
  144. end
  145. if(isempty(A))
  146. error('MATLAB:pca:In1Empty',...
  147. 'Input 1 must not be empty.')
  148. end
  149. %
  150. % Retrieve the dimensions of A.
  151. %
  152. [m n] = size(A);
  153. %
  154. % Check the remaining input arguments.
  155. %
  156. if(size(k,1) ~= 1 || size(k,2) ~= 1)
  157. error('MATLAB:pca:In2Not1x1',...
  158. 'Input 2 must be a scalar.')
  159. end
  160. if(size(its,1) ~= 1 || size(its,2) ~= 1)
  161. error('MATLAB:pca:In3Not1x1',...
  162. 'Input 3 must be a scalar.')
  163. end
  164. if(size(l,1) ~= 1 || size(l,2) ~= 1)
  165. error('MATLAB:pca:In4Not1x1',...
  166. 'Input 4 must be a scalar.')
  167. end
  168. if(k <= 0)
  169. error('MATLAB:pca:In2NonPos',...
  170. 'Input 2 must be > 0.')
  171. end
  172. if((k > m) || (k > n))
  173. error('MATLAB:pca:In2TooBig',...
  174. 'Input 2 must be <= the smallest dimension of Input 1.')
  175. end
  176. if(its < 0)
  177. error('MATLAB:pca:In3Neg',...
  178. 'Input 3 must be >= 0.')
  179. end
  180. if(l < k)
  181. error('MATLAB:pca:In4ltIn2',...
  182. 'Input 4 must be >= Input 2.')
  183. end
  184. %
  185. % SVD A directly if (its+1)*l >= m/1.25 or (its+1)*l >= n/1.25.
  186. %
  187. if(((its+1)*l >= m/1.25) || ((its+1)*l >= n/1.25))
  188. if(~issparse(A))
  189. [U,S,V] = svd(A,'econ');
  190. end
  191. if(issparse(A))
  192. [U,S,V] = svd(full(A),'econ');
  193. end
  194. %
  195. % Retain only the leftmost k columns of U, the leftmost k columns of V,
  196. % and the uppermost leftmost k x k block of S.
  197. %
  198. U = U(:,1:k);
  199. V = V(:,1:k);
  200. S = S(1:k,1:k);
  201. return
  202. end
  203. %
  204. % Take the adjoint of the matrix if m < n.
  205. %
  206. flag = 0;
  207. if(m < n)
  208. A = A';
  209. [m n] = size(A);
  210. flag = 1;
  211. end
  212. %
  213. % Form the adjoint Aa of A.
  214. %
  215. if(its > 0)
  216. Aa = A';
  217. end
  218. %
  219. % Apply A to a random matrix, obtaining H.
  220. %
  221. rand('seed',rand('seed'));
  222. if(isreal(A))
  223. H = A*(2*rand(n,l)-ones(n,l));
  224. end
  225. if(~isreal(A))
  226. H = A*( (2*rand(n,l)-ones(n,l)) + i*(2*rand(n,l)-ones(n,l)) );
  227. end
  228. rand('twister',rand('twister'));
  229. %
  230. % Initialize F to its final size and fill its leftmost block with H.
  231. %
  232. F = zeros(m,(its+1)*l);
  233. F(1:m, 1:l) = H;
  234. %
  235. % Apply A*A' to H a total of its times,
  236. % augmenting F with the new H each time.
  237. %
  238. for it = 1:its
  239. H = Aa*H;
  240. H = A*H;
  241. F(1:m, (1+it*l):((it+1)*l)) = H;
  242. end
  243. clear Aa H;
  244. %
  245. % Form a matrix Q whose columns constitute an orthonormal basis
  246. % for the columns of F.
  247. %
  248. [Q,R,E] = qr(F,0);
  249. clear F R E;
  250. %
  251. % SVD Q'*A to obtain approximations to the singular values
  252. % and right singular vectors of A; adjust the left singular vectors
  253. % of Q'*A to approximate the left singular vectors of A.
  254. %
  255. [U2,S,V] = svd(Q'*A,'econ');
  256. U = Q*U2;
  257. clear Q U2;
  258. %
  259. % Retain only the leftmost k columns of U, the leftmost k columns of V,
  260. % and the uppermost leftmost k x k block of S.
  261. %
  262. U = U(:,1:k);
  263. V = V(:,1:k);
  264. S = S(1:k,1:k);
  265. %
  266. % If flag is nonzero, then swap U and V.
  267. %
  268. if(flag ~= 0)
  269. flag = U;
  270. U = V;
  271. V = flag;
  272. end
  273. clear flag;