-
Notifications
You must be signed in to change notification settings - Fork 0
/
randompca.cpp
807 lines (672 loc) · 22.6 KB
/
randompca.cpp
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
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
/*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 3 of the License, or
* (at your option) any later version.
*
* Copyright (C) 2014-2016 Gad Abraham
* All rights reserved.
*/
#include "randompca.h"
#include "util.h"
#include "svdwide.h"
#include "svdtall.h"
MatrixXd make_gaussian(unsigned int rows, unsigned int cols, long seed)
{
boost::random::mt19937 rng;
rng.seed(seed);
boost::random::normal_distribution<double> nrm;
boost::random::variate_generator<boost::random::mt19937&,
boost::random::normal_distribution<double> > randn(rng, nrm);
MatrixXd G(rows, cols);
std::generate(G.data(), G.data() + G.size(), randn);
return G;
}
// Compute median of pairwise distances on sample of size n from the matrix X
// We're sampling with replacement
// Based on http://www.machinedlearnings.com/2013/08/cosplay.html
// Todo: use Boost accummulators for quantiles
// http://boost-sandbox.sourceforge.net/libs/accumulators/doc/html/accumulators/user_s_guide/the_statistical_accumulators_library.html#accumulators.user_s_guide.the_statistical_accumulators_library.p_square_quantile
double median_dist(MatrixXd& X, unsigned int n, long seed, bool verbose)
{
boost::random::mt19937 rng;
rng.seed(seed);
boost::random::uniform_real_distribution<> dist(0, 1);
double prop = (double)X.rows() / n;
verbose && STDOUT << timestamp() <<
" Computing median Euclidean distance (" << n << " samples)" <<
std::endl;
MatrixXd X2(n, X.cols());
if(n < X.rows())
{
verbose && STDOUT << timestamp() << "Sampling" << std::endl;
// Sample n rows from X
for(unsigned int i = 0, k = 0 ; i < X.rows() ; i++)
{
if(dist(rng) < prop)
{
X2.row(k++) = X.row(i);
if(k == n)
break;
}
}
}
else
X2.noalias() = X;
VectorXd norms = X2.array().square().rowwise().sum();
VectorXd ones = VectorXd::Ones(n);
MatrixXd R = norms * ones.transpose();
MatrixXd D = R + R.transpose() - 2 * X2 * X2.transpose();
unsigned int m = D.size();
double *d = D.data();
std::sort(d, d + m);
double med;
if(m % 2 == 0)
med = (d[m / 2 - 1] + d[m / 2]) / 2;
else
med = d[m / 2];
verbose && STDOUT << timestamp() << "Median Euclidean distance: "
<< med << std::endl;
return med;
}
template <typename Derived>
double var(const MatrixBase<Derived>& x)
{
const unsigned int n = x.size();
const double mean = x.mean();
return (x.array() - mean).square().sum() / (n - 1);
}
template <typename DerivedA, typename DerivedB>
RowVectorXd cov(const MatrixBase<DerivedA>& X, const MatrixBase<DerivedB>& Y)
{
const unsigned int n = X.rows();
const RowVectorXd xmean = X.colwise().mean();
const RowVectorXd ymean = Y.colwise().mean();
return (X.rowwise() - xmean).transpose() * (Y.rowwise() - ymean) / (n - 1);
}
ArrayXXd wilks(const ArrayXd& r2, unsigned int n, unsigned int k)
{
ArrayXd lambda = 1 - r2;
boost::math::fisher_f pf(k, n - k - 1);
ArrayXd pval(r2.size());
ArrayXXd res(r2.size(), 3);
res.col(0) = r2.sqrt();
res.col(1) = (1 - lambda) / lambda * (n - k - 1) / k;
for(unsigned int i = 0 ; i < lambda.size() ; i++)
{
//double fstat = (1 - lambda(i)) / lambda(i) * (n - k - 1) / k;
res(i, 2) = cdf(complement(pf, res(i, 1)));
}
return res;
}
void RandomPCA::pca_fast(MatrixXd& X, unsigned int block_size,
unsigned int ndim, unsigned int maxiter,
double tol, long seed, bool do_loadings)
{
unsigned int N, p;
X_meansd = standardise(X, stand_method_x, verbose);
N = X.rows();
p = X.cols();
SVDWide op(X, verbose);
Spectra::SymEigsSolver<double,
Spectra::LARGEST_ALGE, SVDWide> eigs(&op, ndim, ndim * 2 + 1);
eigs.init();
eigs.compute(maxiter, tol);
double div = 1;
if(divisor == DIVISOR_N1)
div = N - 1;
else if(divisor == DIVISOR_P)
div = p;
if(eigs.info() == Spectra::SUCCESSFUL)
{
U = eigs.eigenvectors();
// Note: _eigenvalues_, not singular values
d = eigs.eigenvalues().array() / div;
if(do_loadings)
{
VectorXd s = d.array().sqrt().inverse() / sqrt(div);
V.noalias() = X.transpose() * U * s.asDiagonal();
}
trace = X.array().square().sum() / div;
pve = d / trace;
Px = U * d.array().sqrt().matrix().asDiagonal();
verbose && STDOUT << timestamp() << "GRM trace: " << trace << std::endl;
}
else
{
throw new std::runtime_error(
std::string("Spectra eigen-decomposition was not successful")
+ ", status: " + std::to_string(eigs.info()));
}
}
void RandomPCA::pca_fast(Data& dat, unsigned int block_size,
unsigned int ndim, unsigned int maxiter, double tol,
long seed, bool do_loadings)
{
unsigned int N = dat.N, p = dat.nsnps;
SVDWideOnline op(dat, block_size, stand_method_x, verbose);
Spectra::SymEigsSolver<double, Spectra::LARGEST_ALGE,
SVDWideOnline> eigs(&op, ndim, ndim * 2 + 1);
eigs.init();
eigs.compute(maxiter, tol);
double div = 1;
if(divisor == DIVISOR_N1)
div = N - 1;
else if(divisor == DIVISOR_P)
div = p;
if(eigs.info() == Spectra::SUCCESSFUL)
{
U = eigs.eigenvectors();
// Note: _eigenvalues_, not singular values
d = eigs.eigenvalues().array() / div;
if(do_loadings)
{
V = MatrixXd::Zero(dat.nsnps, U.cols());
verbose && STDOUT << "Computing loadings" << std::endl;
VectorXd v(dat.nsnps);
for(unsigned int j = 0 ; j < U.cols() ; j++)
{
verbose && STDOUT << "loading " << j << std::endl;
VectorXd u = U.col(j);
op.crossprod(u.data(), v.data());
double s = d(j);
V.col(j) = v * (1.0 / sqrt(s)) / sqrt(div);
}
}
trace = op.trace / div;
pve = d / trace;
Px = U * d.array().sqrt().matrix().asDiagonal();
X_meansd = dat.X_meansd; // TODO: duplication
verbose && STDOUT << timestamp() << "GRM trace: " << trace << std::endl;
}
else
{
throw new std::runtime_error(
std::string("Spectra eigen-decomposition was not successful")
+ ", status: " + std::to_string(eigs.info()));
}
}
double inline sign_scalar(double x)
{
return (0 < x) - (x < 0);
}
VectorXd inline soft_thresh(VectorXd& a, double b)
{
VectorXd s = a.unaryExpr(std::ptr_fun(sign_scalar));
VectorXd d = a.array().abs() - b;
VectorXd z = (d.array() < 0).select(0, d);
return s.array() * z.array();
}
VectorXd norm_thresh(VectorXd& x, double lambda)
{
double s = x.norm();
if(s > 0)
{
x = x.array() / s;
x = soft_thresh(x, lambda);
s = x.norm();
if(s > 0)
x = x.array() / s;
}
return x;
}
void scca_lowmem(MatrixXd& X, MatrixXd &Y, MatrixXd& U, MatrixXd& V,
VectorXd& d, double lambda1, double lambda2,
unsigned int maxiter, double tol, bool verbose)
{
verbose && STDOUT << timestamp()
<< "[scca_lowmem] " << std::endl;
// TODO: X2 and Y2 take up lots of memory
MatrixXd X2 = MatrixXd::Zero(X.rows() + U.cols(), X.cols());
MatrixXd Y2 = MatrixXd::Zero(Y.rows() + U.cols(), Y.cols());
X2.block(0, 0, X.rows(), X.cols()) = X;
Y2.block(0, 0, Y.rows(), Y.cols()) = Y;
VectorXd u, v, u_old, v_old;
for(unsigned int j = 0 ; j < U.cols() ; j++)
{
if(j > 0)
{
X2.row(X.rows() + j) = sqrt(d[j - 1]) * U.col(j - 1).transpose();
Y2.row(Y.rows() + j) = -sqrt(d[j - 1]) * V.col(j - 1).transpose();
}
unsigned int iter = 0;
for( ; iter < maxiter ; iter++)
{
u_old = u = U.col(j);
v_old = v = V.col(j);
u = X2.transpose() * (Y2 * v);
u = norm_thresh(u, lambda1);
U.col(j) = u;
v = Y2.transpose() * (X2 * U.col(j));
v = norm_thresh(v, lambda2);
V.col(j) = v;
if(iter > 0
&& (v_old.array() - v.array()).abs().maxCoeff() < tol
&& (u_old.array() - u.array()).abs().maxCoeff() < tol)
{
verbose && STDOUT << timestamp() << "dim " << j << " finished in "
<< iter << " iterations" << std::endl;
break;
}
}
if(iter >= maxiter)
{
verbose && STDOUT << timestamp()
<< " SCCA did not converge in " << maxiter << " iterations" <<
std::endl;
}
long long nzu = (U.col(j).array() != 0).count();
long long nzv = (V.col(j).array() != 0).count();
verbose && STDOUT << timestamp() << "U_" << j
<< " non-zeros: " << nzu << ", V_" << j
<< " non-zeros: " << nzv << std::endl;
// Use X and Y, not X2 and Y2
d[j] = (X * U.col(j)).transpose() * (Y * V.col(j));
verbose && STDOUT << timestamp() << "d[" << j << "]: "
<< d[j] << std::endl;
}
}
void scca_highmem(MatrixXd& X, MatrixXd &Y, MatrixXd& U, MatrixXd& V,
VectorXd& d, double lambda1, double lambda2,
unsigned int maxiter, double tol, bool verbose)
{
verbose && STDOUT << timestamp()
<< "[scca_highmem] Begin computing X^T Y" << std::endl;
MatrixXd XY = X.transpose() * Y;
verbose && STDOUT << timestamp()
<< "[scaa_highmem] End computing X^T Y" << std::endl;
MatrixXd XYj;
VectorXd u, v, u_old, v_old;
for(unsigned int j = 0 ; j < U.cols() ; j++)
{
verbose && STDOUT << timestamp() << "dim " << j << std::endl;
if(j == 0)
XYj = XY;
else
XYj = XYj - d[j - 1] * U.col(j - 1) * V.col(j - 1).transpose();
unsigned int iter = 0;
for(; iter < maxiter ; iter++)
{
u_old = u = U.col(j);
v_old = v = V.col(j);
u = XYj * v;
u = norm_thresh(u, lambda1);
U.col(j) = u;
v = XYj.transpose() * U.col(j);
v = norm_thresh(v, lambda2);
V.col(j) = v;
if(iter > 0
&& (v_old.array() - v.array()).abs().maxCoeff() < tol
&& (u_old.array() - u.array()).abs().maxCoeff() < tol)
{
verbose && STDOUT << timestamp() << "dim " << j << " finished in "
<< iter << " iterations" << std::endl;
break;
}
}
if(iter >= maxiter)
{
verbose && STDOUT << timestamp()
<< "SCCA did not converge in " << maxiter << " iterations" <<
std::endl;
}
long long nzu = (U.col(j).array() != 0).count();
long long nzv = (V.col(j).array() != 0).count();
verbose && STDOUT << timestamp() << "U_" << j
<< " non-zeros: " << nzu << ", V_" << j
<< " non-zeros: " << nzv << std::endl;
d[j] = U.col(j).transpose() * XYj * V.col(j);
}
}
void RandomPCA::scca(MatrixXd &X, MatrixXd &Y, double lambda1, double lambda2,
long seed, unsigned int ndim, int mem, unsigned int maxiter, double tol)
{
unsigned int k = Y.cols();
MatrixXd M = make_gaussian(k, ndim, seed);
this->scca(X, Y, lambda1, lambda2, seed, ndim, mem, maxiter, tol, M);
}
void RandomPCA::scca(MatrixXd &X, MatrixXd &Y, double lambda1, double lambda2,
long seed, unsigned int ndim, int mem, unsigned int maxiter, double tol,
MatrixXd &V0)
{
X_meansd = standardise(X, stand_method_x);
Y_meansd = standardise(Y, stand_method_y);
verbose && STDOUT << timestamp() << "dim(X): " << dim(X) << std::endl;
verbose && STDOUT << timestamp() << "dim(Y): " << dim(Y) << std::endl;
verbose && STDOUT << timestamp() << "lambda1: " << lambda1
<< " lambda2: " << lambda2 << std::endl;
unsigned int p = X.cols();
//V = make_gaussian(k, ndim, seed);
this->V0 = V0;
V = V0;
U = MatrixXd::Zero(p, ndim);
d = VectorXd::Zero(ndim);
if(mem == HIGHMEM)
scca_highmem(X, Y, U, V, d, lambda1, lambda2, maxiter, tol, verbose);
else
scca_lowmem(X, Y, U, V, d, lambda1, lambda2, maxiter, tol, verbose);
Px = X * U;
Py = Y * V;
}
void RandomPCA::scca(Data &dat, double lambda1, double lambda2,
long seed, unsigned int ndim, int mem, unsigned int maxiter, double tol,
unsigned int block_size)
{
unsigned int k = dat.Y.cols();
MatrixXd M = make_gaussian(k, ndim, seed);
this->scca(dat, lambda1, lambda2, seed, ndim, mem, maxiter, tol, block_size, M);
}
// Block loading of X (genotypes)
// Assumes that data.Y has been set
void RandomPCA::scca(Data &dat, double lambda1, double lambda2,
long seed, unsigned int ndim, int mem, unsigned int maxiter, double tol,
unsigned int block_size, MatrixXd &V0)
{
Y_meansd = standardise(dat.Y, stand_method_y);
verbose && STDOUT << timestamp() << "dim(Y): " << dim(dat.Y) << std::endl;
verbose && STDOUT << timestamp() << "lambda1: " << lambda1
<< " lambda2: " << lambda2 << std::endl;
SVDWideOnline op(dat, block_size, stand_method_x, verbose);
unsigned int p = dat.nsnps;
this->V0 = V0;
V = V0;
U = MatrixXd::Zero(p, ndim);
d = VectorXd::Zero(ndim);
VectorXd uj, vj, uj_old, vj_old;
for(unsigned int j = 0 ; j < U.cols() ; j++)
{
unsigned int iter = 0;
for( ; iter < maxiter ; iter++)
{
uj_old = uj = U.col(j);
vj_old = vj = V.col(j);
// u = X.transpose() * (Y * v);
MatrixXd Yvj = dat.Y * vj;
uj = op.crossprod2(Yvj);
// deflate u
if(j > 0)
{
uj -= U.leftCols(j) * d.head(j - 1).asDiagonal()
* V.leftCols(j).transpose() * vj;
}
uj = norm_thresh(uj, lambda1);
U.col(j) = uj;
// v = Y.transpose() * (X * u);
VectorXd Xuj = op.prod3(uj);
vj = dat.Y.transpose() * Xuj;
// deflate v
if(j > 0)
{
vj -= V.leftCols(j) * d.head(j - 1).asDiagonal()
* U.leftCols(j).transpose() * uj;
}
vj = norm_thresh(vj, lambda2);
V.col(j) = vj;
if(iter > 0
&& (vj_old.array() - vj.array()).abs().maxCoeff() < tol
&& (uj_old.array() - uj.array()).abs().maxCoeff() < tol)
{
verbose && STDOUT << timestamp() << "dim "
<< j << " finished in "
<< iter << " iterations" << std::endl;
break;
}
}
if(iter >= maxiter)
{
verbose && STDOUT << timestamp()
<< " SCCA did not converge in " << maxiter
<< " iterations" << std::endl;
}
long long nzu = (U.col(j).array() != 0).count();
long long nzv = (V.col(j).array() != 0).count();
verbose && STDOUT << timestamp() << "U_" << j
<< " non-zeros: " << nzu << ", V_" << j
<< " non-zeros: " << nzv << std::endl;
VectorXd Xuj = op.prod3(U.col(j));
d[j] = Xuj.transpose() * (dat.Y * V.col(j));
verbose && STDOUT << timestamp() << "d[" << j << "]: "
<< d[j] << std::endl;
}
Px = op.prod3(U);
Py = dat.Y * V;
}
// Single-SNP CCA (like plink.multivariate), offline version (loading all SNPs
// into memory)
void RandomPCA::ucca(MatrixXd &X, MatrixXd &Y)
{
X_meansd = standardise(X, stand_method_x);
Y_meansd = standardise(Y, stand_method_y);
unsigned int n = X.rows();
unsigned int p = X.cols();
double varx;
RowVectorXd covXY;
JacobiSVD<MatrixXd> svd(Y, ComputeThinU | ComputeThinV);
ArrayXd d = svd.singularValues();
MatrixXd V = svd.matrixV();
ArrayXd s(V.cols());
ArrayXd r2 = ArrayXd(p);
for(unsigned int j = 0 ; j < p ; j++)
{
varx = var(X.col(j));
covXY = cov(X.col(j), Y);
// take absolute value to prevent numerical issues with negative numbers
// close to zero
// r2(j) = fabs(1.0 / varx * covXY * covYinv * covXY.transpose());
// covXY * covYinv * covXY'
// = covXY * (V * D^(-2) * V') * covXY'
// = ss'
// where s = covXY * V * D^(-1)
//
// The factor of sqrt(n - 1) comes from the fact that we did SVD of Y,
// but cov(Y) is 1/(n-1) Y'Y when Y is standardised.
s = covXY * V * std::sqrt(double(n - 1));
r2(j) = std::abs((s / d.array()).square().sum() / varx);
}
res = wilks(r2, X.rows(), Y.cols());
}
// Single-SNP CCA (like plink.multivariate), online version (loading one SNP
// at a time)
//
// Assumes data.Y has been set
void RandomPCA::ucca(Data& data)
{
Y_meansd = standardise(data.Y, stand_method_y);
unsigned int n = data.N;
unsigned int p = data.nsnps;
double varx;
RowVectorXd covXY;
verbose && STDOUT << timestamp()
<< "UCCA online mode, N=" << n << " p=" << p << std::endl;
// QR might be faster here
JacobiSVD<MatrixXd> svd(data.Y, ComputeThinU | ComputeThinV);
ArrayXd d = svd.singularValues();
MatrixXd V = svd.matrixV();
ArrayXd s(V.cols());
ArrayXd r2 = ArrayXd(p);
for(unsigned int j = 0 ; j < p ; j++)
{
// No need to explicitly standardise X, since read_snp_block will
// already standardise it internally, assuming that
// data.stand_method_x has been set previously
data.read_snp_block(j, j, false, true);
varx = var(data.X.col(0));
covXY = cov(data.X.col(0), data.Y);
// take absolute value to prevent numerical issues with negative numbers
// that are very close to zero
// r2(j) = fabs(1.0 / varx * covXY * covYinv * covXY.transpose());
// covXY * covYinv * covXY'
// = covXY * (V * D^(-2) * V') * covXY'
// = ss'
// where s = covXY * V * D^(-1)
//
// The factor of sqrt(n - 1) comes from the fact that we did SVD of Y,
// but cov(Y) is 1/(n-1) Y'Y when Y is standardised.
s = covXY * V * std::sqrt(double(n - 1));
r2(j) = std::abs((s / d.array()).square().sum() / varx);
}
res = wilks(r2, n, data.Y.cols());
}
void RandomPCA::check(Data& dat, unsigned int block_size,
std::string evec_file, std::string eval_file)
{
// Read eigenvalues
// Expects no header, no rownames, one eigenvalue per row
verbose && STDOUT << timestamp() << "Loading eigenvalue file '"
<< eval_file << "'" << std::endl;
NamedMatrixWrapper M1 = read_text(eval_file.c_str(), 1, -1, 0);
MatrixXd ev = M1.X;
if(ev.rows() == 0)
throw std::runtime_error("No eigenvalues found in file");
VectorXd eval = ev.col(0);
// Read eigenvectors
// Expects header (colnames), FID and IID cols
verbose && STDOUT << timestamp() << "Loading eigenvector file '"
<< evec_file << "'" << std::endl;
NamedMatrixWrapper M2 = read_text(evec_file.c_str(), 3, -1, 1);
MatrixXd& evec = M2.X;
if(evec.rows() != dat.N)
throw std::runtime_error(
std::string("Eigenvector dimension doesn't match data dimension")
+ " (evec.rows = " + std::to_string(evec.rows())
+ "; dat.N = " + std::to_string(dat.N) + ")");
if(eval.size() != evec.cols())
throw std::runtime_error(
"Eigenvector dimension doesn't match the number of eigenvalues");
check(dat, block_size, evec, eval);
}
// Check the eigenvalues/eigenvectors, computing the root mean squared error
// of E = X X' U / div - U D^2, i.e., averaged over the n * k dimensions.
// assumes WIDE
void RandomPCA::check(Data& dat, unsigned int block_size,
MatrixXd& evec, VectorXd& eval)
{
SVDWideOnline op(dat, block_size, 1, verbose);
unsigned int K = std::min(evec.cols(), eval.size());
// X X' U / div = U D^2
verbose && STDOUT << timestamp()
<< "Checking mean square error between (X X' U) / div and (U D^2)"
<< " for " << K << " dimensions"
<< std::endl;
double div = 1;
if(divisor == DIVISOR_N1)
div = dat.N - 1;
else if(divisor == DIVISOR_P)
div = dat.nsnps;
MatrixXd XXU = op.perform_op_mat(evec);
XXU /= div;
MatrixXd UD2 = evec * eval.asDiagonal();
RowVectorXd rerr = (XXU - UD2).colwise().squaredNorm();
err = rerr.transpose();
for(unsigned int j = 0 ; j < K ; j++)
{
verbose && STDOUT << timestamp() << "eval(" << (j + 1)
<< "): " << eval(j) << ", sum squared error: "
<< err(j) << std::endl;
}
mse = err.sum() / (dat.N * K);
rmse = std::sqrt(mse);
verbose && STDOUT << timestamp() << "Mean squared error: " << mse
<< ", Root mean squared error: " << rmse
<< " (n=" << dat.N << ")" << std::endl;
}
void RandomPCA::check(MatrixXd& X, MatrixXd& evec, VectorXd& eval)
{
X_meansd = standardise(X, stand_method_x, verbose);
unsigned int K = std::min(evec.cols(), eval.size());
unsigned int N = X.rows();
unsigned int p = X.cols();
// X X' U / div = U D^2
verbose && STDOUT << timestamp()
<< "Checking mean square error between (X X' U) and (U D^2)"
<< " for " << K << " dimensions"
<< std::endl;
double div = 1;
if(divisor == DIVISOR_N1)
div = N - 1;
else if(divisor == DIVISOR_P)
div = p;
MatrixXd XXU = X * (X.transpose() * evec) / div;
MatrixXd UD2 = evec * eval.asDiagonal();
RowVectorXd rerr = (XXU - UD2).colwise().squaredNorm();
err = rerr.transpose();
for(unsigned int j = 0 ; j < K ; j++)
{
verbose && STDOUT << timestamp() << "eval(" << (j + 1)
<< "): " << eval(j) << ", sum squared error: "
<< err(j) << std::endl;
}
mse = err.sum() / (N * K);
rmse = std::sqrt(mse);
verbose && STDOUT << timestamp() << "Mean squared error: " << mse
<< ", Root mean squared error: " << rmse
<< " (n=" << N << ")" << std::endl;
}
MatrixXd maf2meansd(MatrixXd maf)
{
MatrixXd X_meansd(maf.rows(), 2);
X_meansd.col(0) = maf * 2.0;
X_meansd.col(1) = maf.array() * 2.0 * (1.0 - maf.array());
return X_meansd;
}
void RandomPCA::project(Data& dat, unsigned int block_size,
std::string loadings_file, std::string maf_file,
std::string meansd_file)
{
// Read the loadings
// TODO: check that SNP ids match
NamedMatrixWrapper M = read_text(loadings_file.c_str(), 3, -1, 1);
V = M.X;
// Read the means+sds or the MAF (and convert MAF to means+sds)
if(maf_file != "")
{
// TODO: missing/non-numeric values?
verbose && STDOUT << timestamp() << "Reading MAF file "
<< maf_file << std::endl;
NamedMatrixWrapper M2 = read_text(maf_file.c_str(), 3, -1, 1);
dat.X_meansd = maf2meansd(M2.X);
dat.use_preloaded_maf = true;
}
else if(meansd_file != "")
{
verbose && STDOUT << timestamp()
<< " Reading mean/stdev file " << meansd_file << std::endl;
NamedMatrixWrapper M2 = read_text(meansd_file.c_str(), 3, -1, 1);
dat.X_meansd = M2.X;
dat.use_preloaded_maf = true;
}
else
{
verbose && STDOUT << timestamp()
<< " Using MAF from the data" << std::endl;
dat.use_preloaded_maf = false;
}
project(dat, block_size);
}
// Project new samples onto existing principal components.
//
// Doesn't do a lot of sanity checking.
//
// Assumes:
// - The loadings matrix V must have been set already
// - dat.X_meansd has beeen set
// - dat.use_preloaded_maf has been set, if needed
void RandomPCA::project(Data& dat, unsigned int block_size)
{
// Check that the SNPs in the data match
SVDWideOnline op(dat, block_size, 1, verbose);
unsigned int k = V.cols();
Px = MatrixXd::Zero(dat.N, k);
double div = 1;
if(divisor == DIVISOR_N1)
div = dat.N - 1;
else if(divisor == DIVISOR_P)
div = V.rows();
VectorXd pxi(dat.N);
for(unsigned int i = 0 ; i < k ; i++)
{
MatrixXd v = V.col(i);
op.prod(v.data(), pxi.data());
Px.col(i) = pxi.array() / sqrt(div); // X V = U D
}
}