irlba
A C++ library for IRLBA
Loading...
Searching...
No Matches
pca.hpp
Go to the documentation of this file.
1#ifndef IRLBA_PCA_HPP
2#define IRLBA_PCA_HPP
3
4#include <cmath>
5#include <utility>
6#include <stdexcept>
7
8#include "compute.hpp"
9#include "Matrix/simple.hpp"
10#include "Matrix/centered.hpp"
11#include "Matrix/scaled.hpp"
12
13#include "Eigen/Dense"
14
20namespace irlba {
21
55template<class InputEigenMatrix_, class OutputEigenMatrix_, class EigenVector_>
56std::pair<bool, int> pca(
57 const InputEigenMatrix_& matrix,
58 bool center,
59 bool scale,
60 Eigen::Index number,
61 OutputEigenMatrix_& scores,
62 OutputEigenMatrix_& rotation,
63 EigenVector_& variances,
64 const Options& options
65) {
66 const Eigen::Index nr = matrix.rows();
67 const Eigen::Index nc = matrix.cols();
68
69 EigenVector_ center0;
70 if (center) {
71 if (nr < 1) {
72 throw std::runtime_error("cannot center with no observations");
73 }
74 center0.resize(nc);
75 }
76
77 EigenVector_ scale0;
78 if (scale) {
79 if (nr < 2) {
80 throw std::runtime_error("cannot scale with fewer than two observations");
81 }
82 scale0.resize(nc);
83 }
84
85 if (center || scale) {
86 EigenVector_ buffer;
87 for (Eigen::Index i = 0; i < nc; ++i) {
88 buffer = matrix.col(i);
89
90 typename EigenVector_::Scalar mean = 0;
91 if (center) {
92 mean = buffer.sum() / nr;
93 center0[i] = mean;
94 }
95
96 if (scale) {
97 typename EigenVector_::Scalar var = 0;
98 for (const auto x : buffer) {
99 var += (x - mean) * (x - mean);
100 }
101
102 if (var) {
103 scale0[i] = std::sqrt(var/(nr - 1));
104 } else {
105 scale0[i] = 1;
106 }
107 }
108 }
109 }
110
111 // Reduce the size of the binary by forcing all compute() calls to use virtual dispatch,
112 // rather than realizing four new instances of the same function with different Matrix subclasses.
113 std::unique_ptr<Matrix<EigenVector_, OutputEigenMatrix_> > ptr;
114 ptr.reset(new SimpleMatrix<EigenVector_, OutputEigenMatrix_, I<decltype(&matrix)> >(&matrix));
115
116 if (center) {
117 std::unique_ptr<Matrix<EigenVector_, OutputEigenMatrix_> > alt;
118 alt.reset(new CenteredMatrix<EigenVector_, OutputEigenMatrix_, I<decltype(ptr)>, I<decltype(&center0)> >(std::move(ptr), &center0));
119 ptr.swap(alt);
120 }
121
122 if (scale) {
123 std::unique_ptr<Matrix<EigenVector_, OutputEigenMatrix_> > alt;
124 alt.reset(new ScaledMatrix<EigenVector_, OutputEigenMatrix_, I<decltype(ptr)>, I<decltype(&scale0)> >(std::move(ptr), &scale0, true, true));
125 ptr.swap(alt);
126 }
127
128 const auto stats = compute(*ptr, number, scores, rotation, variances, options);
129 if (nc > 1) {
130 const auto denom = nc - 1;
131 for (auto& v : variances) {
132 v = v * v / denom;
133 }
134 }
135
136 return stats;
137}
138
144template<class EigenMatrix_, class EigenVector_>
152 EigenMatrix_ scores;
153
160 EigenMatrix_ rotation;
161
166 EigenVector_ variances;
167
172
177};
178
195template<class OutputEigenMatrix_ = Eigen::MatrixXd, class EigenVector_ = Eigen::VectorXd, class InputEigenMatrix_>
196PcaResults<OutputEigenMatrix_, EigenVector_> pca(const InputEigenMatrix_& matrix, bool center, bool scale, Eigen::Index number, const Options& options) {
198 const auto stats = pca(matrix, center, scale, number, output.scores, output.rotation, output.variances, options);
199 output.converged = stats.first;
200 output.iterations = stats.second;
201 return output;
202}
203
204}
205
206#endif
Deferred centering for a matrix.
Deferred centering of a matrix.
Definition centered.hpp:145
Deferred scaling of a matrix.
Definition scaled.hpp:201
A Matrix-compatible wrapper around a "simple" matrix.
Definition simple.hpp:125
Compute an approximate SVD with IRLBA.
Implements IRLBA for approximate SVD.
Definition compute.hpp:22
std::pair< bool, int > compute(const Matrix_ &matrix, const Eigen::Index number, EigenMatrix_ &outU, EigenMatrix_ &outV, EigenVector_ &outD, const Options &options)
Definition compute.hpp:142
std::pair< bool, int > pca(const InputEigenMatrix_ &matrix, bool center, bool scale, Eigen::Index number, OutputEigenMatrix_ &scores, OutputEigenMatrix_ &rotation, EigenVector_ &variances, const Options &options)
Definition pca.hpp:56
Deferred scaling for a matrix.
Simple wrapper around an Eigen matrix.
Options for the IRLBA algorithm.
Definition Options.hpp:16
Results of the IRLBA-based PCA by pca().
Definition pca.hpp:145
EigenMatrix_ scores
Definition pca.hpp:152
EigenVector_ variances
Definition pca.hpp:166
EigenMatrix_ rotation
Definition pca.hpp:160
int iterations
Definition pca.hpp:171
bool converged
Definition pca.hpp:176