Chris@49
|
1 // Copyright (C) 2010-2012 NICTA (www.nicta.com.au)
|
Chris@49
|
2 // Copyright (C) 2010-2012 Conrad Sanderson
|
Chris@49
|
3 // Copyright (C) 2010 Dimitrios Bouzas
|
Chris@49
|
4 // Copyright (C) 2011 Stanislav Funiak
|
Chris@49
|
5 //
|
Chris@49
|
6 // This Source Code Form is subject to the terms of the Mozilla Public
|
Chris@49
|
7 // License, v. 2.0. If a copy of the MPL was not distributed with this
|
Chris@49
|
8 // file, You can obtain one at http://mozilla.org/MPL/2.0/.
|
Chris@49
|
9
|
Chris@49
|
10
|
Chris@49
|
11 //! \addtogroup op_princomp
|
Chris@49
|
12 //! @{
|
Chris@49
|
13
|
Chris@49
|
14
|
Chris@49
|
15
|
Chris@49
|
16 //! \brief
|
Chris@49
|
17 //! principal component analysis -- 4 arguments version
|
Chris@49
|
18 //! computation is done via singular value decomposition
|
Chris@49
|
19 //! coeff_out -> principal component coefficients
|
Chris@49
|
20 //! score_out -> projected samples
|
Chris@49
|
21 //! latent_out -> eigenvalues of principal vectors
|
Chris@49
|
22 //! tsquared_out -> Hotelling's T^2 statistic
|
Chris@49
|
23 template<typename T1>
|
Chris@49
|
24 inline
|
Chris@49
|
25 bool
|
Chris@49
|
26 op_princomp::direct_princomp
|
Chris@49
|
27 (
|
Chris@49
|
28 Mat<typename T1::elem_type>& coeff_out,
|
Chris@49
|
29 Mat<typename T1::elem_type>& score_out,
|
Chris@49
|
30 Col<typename T1::elem_type>& latent_out,
|
Chris@49
|
31 Col<typename T1::elem_type>& tsquared_out,
|
Chris@49
|
32 const Base<typename T1::elem_type, T1>& X,
|
Chris@49
|
33 const typename arma_not_cx<typename T1::elem_type>::result* junk
|
Chris@49
|
34 )
|
Chris@49
|
35 {
|
Chris@49
|
36 arma_extra_debug_sigprint();
|
Chris@49
|
37 arma_ignore(junk);
|
Chris@49
|
38
|
Chris@49
|
39 typedef typename T1::elem_type eT;
|
Chris@49
|
40
|
Chris@49
|
41 const unwrap_check<T1> Y( X.get_ref(), score_out );
|
Chris@49
|
42 const Mat<eT>& in = Y.M;
|
Chris@49
|
43
|
Chris@49
|
44 const uword n_rows = in.n_rows;
|
Chris@49
|
45 const uword n_cols = in.n_cols;
|
Chris@49
|
46
|
Chris@49
|
47 if(n_rows > 1) // more than one sample
|
Chris@49
|
48 {
|
Chris@49
|
49 // subtract the mean - use score_out as temporary matrix
|
Chris@49
|
50 score_out = in - repmat(mean(in), n_rows, 1);
|
Chris@49
|
51
|
Chris@49
|
52 // singular value decomposition
|
Chris@49
|
53 Mat<eT> U;
|
Chris@49
|
54 Col<eT> s;
|
Chris@49
|
55
|
Chris@49
|
56 const bool svd_ok = svd(U,s,coeff_out,score_out);
|
Chris@49
|
57
|
Chris@49
|
58 if(svd_ok == false)
|
Chris@49
|
59 {
|
Chris@49
|
60 return false;
|
Chris@49
|
61 }
|
Chris@49
|
62
|
Chris@49
|
63
|
Chris@49
|
64 //U.reset(); // TODO: do we need this ? U will get automatically deleted anyway
|
Chris@49
|
65
|
Chris@49
|
66 // normalize the eigenvalues
|
Chris@49
|
67 s /= std::sqrt( double(n_rows - 1) );
|
Chris@49
|
68
|
Chris@49
|
69 // project the samples to the principals
|
Chris@49
|
70 score_out *= coeff_out;
|
Chris@49
|
71
|
Chris@49
|
72 if(n_rows <= n_cols) // number of samples is less than their dimensionality
|
Chris@49
|
73 {
|
Chris@49
|
74 score_out.cols(n_rows-1,n_cols-1).zeros();
|
Chris@49
|
75
|
Chris@49
|
76 //Col<eT> s_tmp = zeros< Col<eT> >(n_cols);
|
Chris@49
|
77 Col<eT> s_tmp(n_cols);
|
Chris@49
|
78 s_tmp.zeros();
|
Chris@49
|
79
|
Chris@49
|
80 s_tmp.rows(0,n_rows-2) = s.rows(0,n_rows-2);
|
Chris@49
|
81 s = s_tmp;
|
Chris@49
|
82
|
Chris@49
|
83 // compute the Hotelling's T-squared
|
Chris@49
|
84 s_tmp.rows(0,n_rows-2) = eT(1) / s_tmp.rows(0,n_rows-2);
|
Chris@49
|
85
|
Chris@49
|
86 const Mat<eT> S = score_out * diagmat(Col<eT>(s_tmp));
|
Chris@49
|
87 tsquared_out = sum(S%S,1);
|
Chris@49
|
88 }
|
Chris@49
|
89 else
|
Chris@49
|
90 {
|
Chris@49
|
91 // compute the Hotelling's T-squared
|
Chris@49
|
92 const Mat<eT> S = score_out * diagmat(Col<eT>( eT(1) / s));
|
Chris@49
|
93 tsquared_out = sum(S%S,1);
|
Chris@49
|
94 }
|
Chris@49
|
95
|
Chris@49
|
96 // compute the eigenvalues of the principal vectors
|
Chris@49
|
97 latent_out = s%s;
|
Chris@49
|
98 }
|
Chris@49
|
99 else // 0 or 1 samples
|
Chris@49
|
100 {
|
Chris@49
|
101 coeff_out.eye(n_cols, n_cols);
|
Chris@49
|
102
|
Chris@49
|
103 score_out.copy_size(in);
|
Chris@49
|
104 score_out.zeros();
|
Chris@49
|
105
|
Chris@49
|
106 latent_out.set_size(n_cols);
|
Chris@49
|
107 latent_out.zeros();
|
Chris@49
|
108
|
Chris@49
|
109 tsquared_out.set_size(n_rows);
|
Chris@49
|
110 tsquared_out.zeros();
|
Chris@49
|
111 }
|
Chris@49
|
112
|
Chris@49
|
113 return true;
|
Chris@49
|
114 }
|
Chris@49
|
115
|
Chris@49
|
116
|
Chris@49
|
117
|
Chris@49
|
118 //! \brief
|
Chris@49
|
119 //! principal component analysis -- 3 arguments version
|
Chris@49
|
120 //! computation is done via singular value decomposition
|
Chris@49
|
121 //! coeff_out -> principal component coefficients
|
Chris@49
|
122 //! score_out -> projected samples
|
Chris@49
|
123 //! latent_out -> eigenvalues of principal vectors
|
Chris@49
|
124 template<typename T1>
|
Chris@49
|
125 inline
|
Chris@49
|
126 bool
|
Chris@49
|
127 op_princomp::direct_princomp
|
Chris@49
|
128 (
|
Chris@49
|
129 Mat<typename T1::elem_type>& coeff_out,
|
Chris@49
|
130 Mat<typename T1::elem_type>& score_out,
|
Chris@49
|
131 Col<typename T1::elem_type>& latent_out,
|
Chris@49
|
132 const Base<typename T1::elem_type, T1>& X,
|
Chris@49
|
133 const typename arma_not_cx<typename T1::elem_type>::result* junk
|
Chris@49
|
134 )
|
Chris@49
|
135 {
|
Chris@49
|
136 arma_extra_debug_sigprint();
|
Chris@49
|
137 arma_ignore(junk);
|
Chris@49
|
138
|
Chris@49
|
139 typedef typename T1::elem_type eT;
|
Chris@49
|
140
|
Chris@49
|
141 const unwrap_check<T1> Y( X.get_ref(), score_out );
|
Chris@49
|
142 const Mat<eT>& in = Y.M;
|
Chris@49
|
143
|
Chris@49
|
144 const uword n_rows = in.n_rows;
|
Chris@49
|
145 const uword n_cols = in.n_cols;
|
Chris@49
|
146
|
Chris@49
|
147 if(n_rows > 1) // more than one sample
|
Chris@49
|
148 {
|
Chris@49
|
149 // subtract the mean - use score_out as temporary matrix
|
Chris@49
|
150 score_out = in - repmat(mean(in), n_rows, 1);
|
Chris@49
|
151
|
Chris@49
|
152 // singular value decomposition
|
Chris@49
|
153 Mat<eT> U;
|
Chris@49
|
154 Col<eT> s;
|
Chris@49
|
155
|
Chris@49
|
156 const bool svd_ok = svd(U,s,coeff_out,score_out);
|
Chris@49
|
157
|
Chris@49
|
158 if(svd_ok == false)
|
Chris@49
|
159 {
|
Chris@49
|
160 return false;
|
Chris@49
|
161 }
|
Chris@49
|
162
|
Chris@49
|
163
|
Chris@49
|
164 // U.reset();
|
Chris@49
|
165
|
Chris@49
|
166 // normalize the eigenvalues
|
Chris@49
|
167 s /= std::sqrt( double(n_rows - 1) );
|
Chris@49
|
168
|
Chris@49
|
169 // project the samples to the principals
|
Chris@49
|
170 score_out *= coeff_out;
|
Chris@49
|
171
|
Chris@49
|
172 if(n_rows <= n_cols) // number of samples is less than their dimensionality
|
Chris@49
|
173 {
|
Chris@49
|
174 score_out.cols(n_rows-1,n_cols-1).zeros();
|
Chris@49
|
175
|
Chris@49
|
176 Col<eT> s_tmp = zeros< Col<eT> >(n_cols);
|
Chris@49
|
177 s_tmp.rows(0,n_rows-2) = s.rows(0,n_rows-2);
|
Chris@49
|
178 s = s_tmp;
|
Chris@49
|
179 }
|
Chris@49
|
180
|
Chris@49
|
181 // compute the eigenvalues of the principal vectors
|
Chris@49
|
182 latent_out = s%s;
|
Chris@49
|
183
|
Chris@49
|
184 }
|
Chris@49
|
185 else // 0 or 1 samples
|
Chris@49
|
186 {
|
Chris@49
|
187 coeff_out.eye(n_cols, n_cols);
|
Chris@49
|
188
|
Chris@49
|
189 score_out.copy_size(in);
|
Chris@49
|
190 score_out.zeros();
|
Chris@49
|
191
|
Chris@49
|
192 latent_out.set_size(n_cols);
|
Chris@49
|
193 latent_out.zeros();
|
Chris@49
|
194 }
|
Chris@49
|
195
|
Chris@49
|
196 return true;
|
Chris@49
|
197 }
|
Chris@49
|
198
|
Chris@49
|
199
|
Chris@49
|
200
|
Chris@49
|
201 //! \brief
|
Chris@49
|
202 //! principal component analysis -- 2 arguments version
|
Chris@49
|
203 //! computation is done via singular value decomposition
|
Chris@49
|
204 //! coeff_out -> principal component coefficients
|
Chris@49
|
205 //! score_out -> projected samples
|
Chris@49
|
206 template<typename T1>
|
Chris@49
|
207 inline
|
Chris@49
|
208 bool
|
Chris@49
|
209 op_princomp::direct_princomp
|
Chris@49
|
210 (
|
Chris@49
|
211 Mat<typename T1::elem_type>& coeff_out,
|
Chris@49
|
212 Mat<typename T1::elem_type>& score_out,
|
Chris@49
|
213 const Base<typename T1::elem_type, T1>& X,
|
Chris@49
|
214 const typename arma_not_cx<typename T1::elem_type>::result* junk
|
Chris@49
|
215 )
|
Chris@49
|
216 {
|
Chris@49
|
217 arma_extra_debug_sigprint();
|
Chris@49
|
218 arma_ignore(junk);
|
Chris@49
|
219
|
Chris@49
|
220 typedef typename T1::elem_type eT;
|
Chris@49
|
221
|
Chris@49
|
222 const unwrap_check<T1> Y( X.get_ref(), score_out );
|
Chris@49
|
223 const Mat<eT>& in = Y.M;
|
Chris@49
|
224
|
Chris@49
|
225 const uword n_rows = in.n_rows;
|
Chris@49
|
226 const uword n_cols = in.n_cols;
|
Chris@49
|
227
|
Chris@49
|
228 if(n_rows > 1) // more than one sample
|
Chris@49
|
229 {
|
Chris@49
|
230 // subtract the mean - use score_out as temporary matrix
|
Chris@49
|
231 score_out = in - repmat(mean(in), n_rows, 1);
|
Chris@49
|
232
|
Chris@49
|
233 // singular value decomposition
|
Chris@49
|
234 Mat<eT> U;
|
Chris@49
|
235 Col<eT> s;
|
Chris@49
|
236
|
Chris@49
|
237 const bool svd_ok = svd(U,s,coeff_out,score_out);
|
Chris@49
|
238
|
Chris@49
|
239 if(svd_ok == false)
|
Chris@49
|
240 {
|
Chris@49
|
241 return false;
|
Chris@49
|
242 }
|
Chris@49
|
243
|
Chris@49
|
244 // U.reset();
|
Chris@49
|
245
|
Chris@49
|
246 // normalize the eigenvalues
|
Chris@49
|
247 s /= std::sqrt( double(n_rows - 1) );
|
Chris@49
|
248
|
Chris@49
|
249 // project the samples to the principals
|
Chris@49
|
250 score_out *= coeff_out;
|
Chris@49
|
251
|
Chris@49
|
252 if(n_rows <= n_cols) // number of samples is less than their dimensionality
|
Chris@49
|
253 {
|
Chris@49
|
254 score_out.cols(n_rows-1,n_cols-1).zeros();
|
Chris@49
|
255
|
Chris@49
|
256 Col<eT> s_tmp = zeros< Col<eT> >(n_cols);
|
Chris@49
|
257 s_tmp.rows(0,n_rows-2) = s.rows(0,n_rows-2);
|
Chris@49
|
258 s = s_tmp;
|
Chris@49
|
259 }
|
Chris@49
|
260 }
|
Chris@49
|
261 else // 0 or 1 samples
|
Chris@49
|
262 {
|
Chris@49
|
263 coeff_out.eye(n_cols, n_cols);
|
Chris@49
|
264 score_out.copy_size(in);
|
Chris@49
|
265 score_out.zeros();
|
Chris@49
|
266 }
|
Chris@49
|
267
|
Chris@49
|
268 return true;
|
Chris@49
|
269 }
|
Chris@49
|
270
|
Chris@49
|
271
|
Chris@49
|
272
|
Chris@49
|
273 //! \brief
|
Chris@49
|
274 //! principal component analysis -- 1 argument version
|
Chris@49
|
275 //! computation is done via singular value decomposition
|
Chris@49
|
276 //! coeff_out -> principal component coefficients
|
Chris@49
|
277 template<typename T1>
|
Chris@49
|
278 inline
|
Chris@49
|
279 bool
|
Chris@49
|
280 op_princomp::direct_princomp
|
Chris@49
|
281 (
|
Chris@49
|
282 Mat<typename T1::elem_type>& coeff_out,
|
Chris@49
|
283 const Base<typename T1::elem_type, T1>& X,
|
Chris@49
|
284 const typename arma_not_cx<typename T1::elem_type>::result* junk
|
Chris@49
|
285 )
|
Chris@49
|
286 {
|
Chris@49
|
287 arma_extra_debug_sigprint();
|
Chris@49
|
288 arma_ignore(junk);
|
Chris@49
|
289
|
Chris@49
|
290 typedef typename T1::elem_type eT;
|
Chris@49
|
291
|
Chris@49
|
292 const unwrap<T1> Y( X.get_ref() );
|
Chris@49
|
293 const Mat<eT>& in = Y.M;
|
Chris@49
|
294
|
Chris@49
|
295 if(in.n_elem != 0)
|
Chris@49
|
296 {
|
Chris@49
|
297 // singular value decomposition
|
Chris@49
|
298 Mat<eT> U;
|
Chris@49
|
299 Col<eT> s;
|
Chris@49
|
300
|
Chris@49
|
301 const Mat<eT> tmp = in - repmat(mean(in), in.n_rows, 1);
|
Chris@49
|
302
|
Chris@49
|
303 const bool svd_ok = svd(U,s,coeff_out, tmp);
|
Chris@49
|
304
|
Chris@49
|
305 if(svd_ok == false)
|
Chris@49
|
306 {
|
Chris@49
|
307 return false;
|
Chris@49
|
308 }
|
Chris@49
|
309 }
|
Chris@49
|
310 else
|
Chris@49
|
311 {
|
Chris@49
|
312 coeff_out.eye(in.n_cols, in.n_cols);
|
Chris@49
|
313 }
|
Chris@49
|
314
|
Chris@49
|
315 return true;
|
Chris@49
|
316 }
|
Chris@49
|
317
|
Chris@49
|
318
|
Chris@49
|
319
|
Chris@49
|
320 //! \brief
|
Chris@49
|
321 //! principal component analysis -- 4 arguments complex version
|
Chris@49
|
322 //! computation is done via singular value decomposition
|
Chris@49
|
323 //! coeff_out -> principal component coefficients
|
Chris@49
|
324 //! score_out -> projected samples
|
Chris@49
|
325 //! latent_out -> eigenvalues of principal vectors
|
Chris@49
|
326 //! tsquared_out -> Hotelling's T^2 statistic
|
Chris@49
|
327 template<typename T1>
|
Chris@49
|
328 inline
|
Chris@49
|
329 bool
|
Chris@49
|
330 op_princomp::direct_princomp
|
Chris@49
|
331 (
|
Chris@49
|
332 Mat< std::complex<typename T1::pod_type> >& coeff_out,
|
Chris@49
|
333 Mat< std::complex<typename T1::pod_type> >& score_out,
|
Chris@49
|
334 Col< typename T1::pod_type >& latent_out,
|
Chris@49
|
335 Col< std::complex<typename T1::pod_type> >& tsquared_out,
|
Chris@49
|
336 const Base< std::complex<typename T1::pod_type>, T1 >& X,
|
Chris@49
|
337 const typename arma_cx_only<typename T1::elem_type>::result* junk
|
Chris@49
|
338 )
|
Chris@49
|
339 {
|
Chris@49
|
340 arma_extra_debug_sigprint();
|
Chris@49
|
341 arma_ignore(junk);
|
Chris@49
|
342
|
Chris@49
|
343 typedef typename T1::pod_type T;
|
Chris@49
|
344 typedef std::complex<T> eT;
|
Chris@49
|
345
|
Chris@49
|
346 const unwrap_check<T1> Y( X.get_ref(), score_out );
|
Chris@49
|
347 const Mat<eT>& in = Y.M;
|
Chris@49
|
348
|
Chris@49
|
349 const uword n_rows = in.n_rows;
|
Chris@49
|
350 const uword n_cols = in.n_cols;
|
Chris@49
|
351
|
Chris@49
|
352 if(n_rows > 1) // more than one sample
|
Chris@49
|
353 {
|
Chris@49
|
354 // subtract the mean - use score_out as temporary matrix
|
Chris@49
|
355 score_out = in - repmat(mean(in), n_rows, 1);
|
Chris@49
|
356
|
Chris@49
|
357 // singular value decomposition
|
Chris@49
|
358 Mat<eT> U;
|
Chris@49
|
359 Col<T> s;
|
Chris@49
|
360
|
Chris@49
|
361 const bool svd_ok = svd(U,s,coeff_out,score_out);
|
Chris@49
|
362
|
Chris@49
|
363 if(svd_ok == false)
|
Chris@49
|
364 {
|
Chris@49
|
365 return false;
|
Chris@49
|
366 }
|
Chris@49
|
367
|
Chris@49
|
368
|
Chris@49
|
369 //U.reset();
|
Chris@49
|
370
|
Chris@49
|
371 // normalize the eigenvalues
|
Chris@49
|
372 s /= std::sqrt( double(n_rows - 1) );
|
Chris@49
|
373
|
Chris@49
|
374 // project the samples to the principals
|
Chris@49
|
375 score_out *= coeff_out;
|
Chris@49
|
376
|
Chris@49
|
377 if(n_rows <= n_cols) // number of samples is less than their dimensionality
|
Chris@49
|
378 {
|
Chris@49
|
379 score_out.cols(n_rows-1,n_cols-1).zeros();
|
Chris@49
|
380
|
Chris@49
|
381 Col<T> s_tmp = zeros< Col<T> >(n_cols);
|
Chris@49
|
382 s_tmp.rows(0,n_rows-2) = s.rows(0,n_rows-2);
|
Chris@49
|
383 s = s_tmp;
|
Chris@49
|
384
|
Chris@49
|
385 // compute the Hotelling's T-squared
|
Chris@49
|
386 s_tmp.rows(0,n_rows-2) = 1.0 / s_tmp.rows(0,n_rows-2);
|
Chris@49
|
387 const Mat<eT> S = score_out * diagmat(Col<T>(s_tmp));
|
Chris@49
|
388 tsquared_out = sum(S%S,1);
|
Chris@49
|
389 }
|
Chris@49
|
390 else
|
Chris@49
|
391 {
|
Chris@49
|
392 // compute the Hotelling's T-squared
|
Chris@49
|
393 const Mat<eT> S = score_out * diagmat(Col<T>(T(1) / s));
|
Chris@49
|
394 tsquared_out = sum(S%S,1);
|
Chris@49
|
395 }
|
Chris@49
|
396
|
Chris@49
|
397 // compute the eigenvalues of the principal vectors
|
Chris@49
|
398 latent_out = s%s;
|
Chris@49
|
399
|
Chris@49
|
400 }
|
Chris@49
|
401 else // 0 or 1 samples
|
Chris@49
|
402 {
|
Chris@49
|
403 coeff_out.eye(n_cols, n_cols);
|
Chris@49
|
404
|
Chris@49
|
405 score_out.copy_size(in);
|
Chris@49
|
406 score_out.zeros();
|
Chris@49
|
407
|
Chris@49
|
408 latent_out.set_size(n_cols);
|
Chris@49
|
409 latent_out.zeros();
|
Chris@49
|
410
|
Chris@49
|
411 tsquared_out.set_size(n_rows);
|
Chris@49
|
412 tsquared_out.zeros();
|
Chris@49
|
413 }
|
Chris@49
|
414
|
Chris@49
|
415 return true;
|
Chris@49
|
416 }
|
Chris@49
|
417
|
Chris@49
|
418
|
Chris@49
|
419
|
Chris@49
|
420 //! \brief
|
Chris@49
|
421 //! principal component analysis -- 3 arguments complex version
|
Chris@49
|
422 //! computation is done via singular value decomposition
|
Chris@49
|
423 //! coeff_out -> principal component coefficients
|
Chris@49
|
424 //! score_out -> projected samples
|
Chris@49
|
425 //! latent_out -> eigenvalues of principal vectors
|
Chris@49
|
426 template<typename T1>
|
Chris@49
|
427 inline
|
Chris@49
|
428 bool
|
Chris@49
|
429 op_princomp::direct_princomp
|
Chris@49
|
430 (
|
Chris@49
|
431 Mat< std::complex<typename T1::pod_type> >& coeff_out,
|
Chris@49
|
432 Mat< std::complex<typename T1::pod_type> >& score_out,
|
Chris@49
|
433 Col< typename T1::pod_type >& latent_out,
|
Chris@49
|
434 const Base< std::complex<typename T1::pod_type>, T1 >& X,
|
Chris@49
|
435 const typename arma_cx_only<typename T1::elem_type>::result* junk
|
Chris@49
|
436 )
|
Chris@49
|
437 {
|
Chris@49
|
438 arma_extra_debug_sigprint();
|
Chris@49
|
439 arma_ignore(junk);
|
Chris@49
|
440
|
Chris@49
|
441 typedef typename T1::pod_type T;
|
Chris@49
|
442 typedef std::complex<T> eT;
|
Chris@49
|
443
|
Chris@49
|
444 const unwrap_check<T1> Y( X.get_ref(), score_out );
|
Chris@49
|
445 const Mat<eT>& in = Y.M;
|
Chris@49
|
446
|
Chris@49
|
447 const uword n_rows = in.n_rows;
|
Chris@49
|
448 const uword n_cols = in.n_cols;
|
Chris@49
|
449
|
Chris@49
|
450 if(n_rows > 1) // more than one sample
|
Chris@49
|
451 {
|
Chris@49
|
452 // subtract the mean - use score_out as temporary matrix
|
Chris@49
|
453 score_out = in - repmat(mean(in), n_rows, 1);
|
Chris@49
|
454
|
Chris@49
|
455 // singular value decomposition
|
Chris@49
|
456 Mat<eT> U;
|
Chris@49
|
457 Col< T> s;
|
Chris@49
|
458
|
Chris@49
|
459 const bool svd_ok = svd(U,s,coeff_out,score_out);
|
Chris@49
|
460
|
Chris@49
|
461 if(svd_ok == false)
|
Chris@49
|
462 {
|
Chris@49
|
463 return false;
|
Chris@49
|
464 }
|
Chris@49
|
465
|
Chris@49
|
466
|
Chris@49
|
467 // U.reset();
|
Chris@49
|
468
|
Chris@49
|
469 // normalize the eigenvalues
|
Chris@49
|
470 s /= std::sqrt( double(n_rows - 1) );
|
Chris@49
|
471
|
Chris@49
|
472 // project the samples to the principals
|
Chris@49
|
473 score_out *= coeff_out;
|
Chris@49
|
474
|
Chris@49
|
475 if(n_rows <= n_cols) // number of samples is less than their dimensionality
|
Chris@49
|
476 {
|
Chris@49
|
477 score_out.cols(n_rows-1,n_cols-1).zeros();
|
Chris@49
|
478
|
Chris@49
|
479 Col<T> s_tmp = zeros< Col<T> >(n_cols);
|
Chris@49
|
480 s_tmp.rows(0,n_rows-2) = s.rows(0,n_rows-2);
|
Chris@49
|
481 s = s_tmp;
|
Chris@49
|
482 }
|
Chris@49
|
483
|
Chris@49
|
484 // compute the eigenvalues of the principal vectors
|
Chris@49
|
485 latent_out = s%s;
|
Chris@49
|
486
|
Chris@49
|
487 }
|
Chris@49
|
488 else // 0 or 1 samples
|
Chris@49
|
489 {
|
Chris@49
|
490 coeff_out.eye(n_cols, n_cols);
|
Chris@49
|
491
|
Chris@49
|
492 score_out.copy_size(in);
|
Chris@49
|
493 score_out.zeros();
|
Chris@49
|
494
|
Chris@49
|
495 latent_out.set_size(n_cols);
|
Chris@49
|
496 latent_out.zeros();
|
Chris@49
|
497 }
|
Chris@49
|
498
|
Chris@49
|
499 return true;
|
Chris@49
|
500 }
|
Chris@49
|
501
|
Chris@49
|
502
|
Chris@49
|
503
|
Chris@49
|
504 //! \brief
|
Chris@49
|
505 //! principal component analysis -- 2 arguments complex version
|
Chris@49
|
506 //! computation is done via singular value decomposition
|
Chris@49
|
507 //! coeff_out -> principal component coefficients
|
Chris@49
|
508 //! score_out -> projected samples
|
Chris@49
|
509 template<typename T1>
|
Chris@49
|
510 inline
|
Chris@49
|
511 bool
|
Chris@49
|
512 op_princomp::direct_princomp
|
Chris@49
|
513 (
|
Chris@49
|
514 Mat< std::complex<typename T1::pod_type> >& coeff_out,
|
Chris@49
|
515 Mat< std::complex<typename T1::pod_type> >& score_out,
|
Chris@49
|
516 const Base< std::complex<typename T1::pod_type>, T1 >& X,
|
Chris@49
|
517 const typename arma_cx_only<typename T1::elem_type>::result* junk
|
Chris@49
|
518 )
|
Chris@49
|
519 {
|
Chris@49
|
520 arma_extra_debug_sigprint();
|
Chris@49
|
521 arma_ignore(junk);
|
Chris@49
|
522
|
Chris@49
|
523 typedef typename T1::pod_type T;
|
Chris@49
|
524 typedef std::complex<T> eT;
|
Chris@49
|
525
|
Chris@49
|
526 const unwrap_check<T1> Y( X.get_ref(), score_out );
|
Chris@49
|
527 const Mat<eT>& in = Y.M;
|
Chris@49
|
528
|
Chris@49
|
529 const uword n_rows = in.n_rows;
|
Chris@49
|
530 const uword n_cols = in.n_cols;
|
Chris@49
|
531
|
Chris@49
|
532 if(n_rows > 1) // more than one sample
|
Chris@49
|
533 {
|
Chris@49
|
534 // subtract the mean - use score_out as temporary matrix
|
Chris@49
|
535 score_out = in - repmat(mean(in), n_rows, 1);
|
Chris@49
|
536
|
Chris@49
|
537 // singular value decomposition
|
Chris@49
|
538 Mat<eT> U;
|
Chris@49
|
539 Col< T> s;
|
Chris@49
|
540
|
Chris@49
|
541 const bool svd_ok = svd(U,s,coeff_out,score_out);
|
Chris@49
|
542
|
Chris@49
|
543 if(svd_ok == false)
|
Chris@49
|
544 {
|
Chris@49
|
545 return false;
|
Chris@49
|
546 }
|
Chris@49
|
547
|
Chris@49
|
548 // U.reset();
|
Chris@49
|
549
|
Chris@49
|
550 // normalize the eigenvalues
|
Chris@49
|
551 s /= std::sqrt( double(n_rows - 1) );
|
Chris@49
|
552
|
Chris@49
|
553 // project the samples to the principals
|
Chris@49
|
554 score_out *= coeff_out;
|
Chris@49
|
555
|
Chris@49
|
556 if(n_rows <= n_cols) // number of samples is less than their dimensionality
|
Chris@49
|
557 {
|
Chris@49
|
558 score_out.cols(n_rows-1,n_cols-1).zeros();
|
Chris@49
|
559 }
|
Chris@49
|
560
|
Chris@49
|
561 }
|
Chris@49
|
562 else // 0 or 1 samples
|
Chris@49
|
563 {
|
Chris@49
|
564 coeff_out.eye(n_cols, n_cols);
|
Chris@49
|
565
|
Chris@49
|
566 score_out.copy_size(in);
|
Chris@49
|
567 score_out.zeros();
|
Chris@49
|
568 }
|
Chris@49
|
569
|
Chris@49
|
570 return true;
|
Chris@49
|
571 }
|
Chris@49
|
572
|
Chris@49
|
573
|
Chris@49
|
574
|
Chris@49
|
575 //! \brief
|
Chris@49
|
576 //! principal component analysis -- 1 argument complex version
|
Chris@49
|
577 //! computation is done via singular value decomposition
|
Chris@49
|
578 //! coeff_out -> principal component coefficients
|
Chris@49
|
579 template<typename T1>
|
Chris@49
|
580 inline
|
Chris@49
|
581 bool
|
Chris@49
|
582 op_princomp::direct_princomp
|
Chris@49
|
583 (
|
Chris@49
|
584 Mat< std::complex<typename T1::pod_type> >& coeff_out,
|
Chris@49
|
585 const Base< std::complex<typename T1::pod_type>, T1 >& X,
|
Chris@49
|
586 const typename arma_cx_only<typename T1::elem_type>::result* junk
|
Chris@49
|
587 )
|
Chris@49
|
588 {
|
Chris@49
|
589 arma_extra_debug_sigprint();
|
Chris@49
|
590 arma_ignore(junk);
|
Chris@49
|
591
|
Chris@49
|
592 typedef typename T1::pod_type T;
|
Chris@49
|
593 typedef std::complex<T> eT;
|
Chris@49
|
594
|
Chris@49
|
595 const unwrap<T1> Y( X.get_ref() );
|
Chris@49
|
596 const Mat<eT>& in = Y.M;
|
Chris@49
|
597
|
Chris@49
|
598 if(in.n_elem != 0)
|
Chris@49
|
599 {
|
Chris@49
|
600 // singular value decomposition
|
Chris@49
|
601 Mat<eT> U;
|
Chris@49
|
602 Col< T> s;
|
Chris@49
|
603
|
Chris@49
|
604 const Mat<eT> tmp = in - repmat(mean(in), in.n_rows, 1);
|
Chris@49
|
605
|
Chris@49
|
606 const bool svd_ok = svd(U,s,coeff_out, tmp);
|
Chris@49
|
607
|
Chris@49
|
608 if(svd_ok == false)
|
Chris@49
|
609 {
|
Chris@49
|
610 return false;
|
Chris@49
|
611 }
|
Chris@49
|
612 }
|
Chris@49
|
613 else
|
Chris@49
|
614 {
|
Chris@49
|
615 coeff_out.eye(in.n_cols, in.n_cols);
|
Chris@49
|
616 }
|
Chris@49
|
617
|
Chris@49
|
618 return true;
|
Chris@49
|
619 }
|
Chris@49
|
620
|
Chris@49
|
621
|
Chris@49
|
622
|
Chris@49
|
623 template<typename T1>
|
Chris@49
|
624 inline
|
Chris@49
|
625 void
|
Chris@49
|
626 op_princomp::apply
|
Chris@49
|
627 (
|
Chris@49
|
628 Mat<typename T1::elem_type>& out,
|
Chris@49
|
629 const Op<T1,op_princomp>& in
|
Chris@49
|
630 )
|
Chris@49
|
631 {
|
Chris@49
|
632 arma_extra_debug_sigprint();
|
Chris@49
|
633
|
Chris@49
|
634 typedef typename T1::elem_type eT;
|
Chris@49
|
635
|
Chris@49
|
636 const unwrap_check<T1> tmp(in.m, out);
|
Chris@49
|
637 const Mat<eT>& A = tmp.M;
|
Chris@49
|
638
|
Chris@49
|
639 const bool status = op_princomp::direct_princomp(out, A);
|
Chris@49
|
640
|
Chris@49
|
641 if(status == false)
|
Chris@49
|
642 {
|
Chris@49
|
643 out.reset();
|
Chris@49
|
644
|
Chris@49
|
645 arma_bad("princomp(): failed to converge");
|
Chris@49
|
646 }
|
Chris@49
|
647 }
|
Chris@49
|
648
|
Chris@49
|
649
|
Chris@49
|
650
|
Chris@49
|
651 //! @}
|