v2.0.0
Loading...
Searching...
No Matches
ica.cpp
Go to the documentation of this file.
1//=============================================================================================================
33
34//=============================================================================================================
35// INCLUDES
36//=============================================================================================================
37
38#include "ica.h"
39
40//=============================================================================================================
41// EIGEN INCLUDES
42//=============================================================================================================
43
44#include <Eigen/Dense>
45
46//=============================================================================================================
47// QT INCLUDES
48//=============================================================================================================
49
50#include <QDebug>
51
52//=============================================================================================================
53// C++ INCLUDES
54//=============================================================================================================
55
56#include <cmath>
57#include <random>
58
59//=============================================================================================================
60// USED NAMESPACES
61//=============================================================================================================
62
63using namespace UTILSLIB;
64using namespace Eigen;
65
66//=============================================================================================================
67// STATIC MEMBER DEFINITIONS
68//=============================================================================================================
69
70//=============================================================================================================
71// PRIVATE HELPERS
72//=============================================================================================================
73
74namespace {
75
76//=============================================================================================================
80void gramSchmidt(VectorXd& w, const MatrixXd& W, int iComp)
81{
82 for (int j = 0; j < iComp; ++j) {
83 w -= w.dot(W.row(j)) * W.row(j).transpose();
84 }
85}
86
87} // anonymous namespace
88
89//=============================================================================================================
90// MEMBER DEFINITIONS
91//=============================================================================================================
92
93IcaResult ICA::run(const MatrixXd& matData,
94 int nComponents,
95 int maxIter,
96 double tol,
97 int randomSeed)
98{
99 const int nCh = static_cast<int>(matData.rows());
100 const int nSamples = static_cast<int>(matData.cols());
101
102 if (nComponents <= 0 || nComponents > nCh) {
103 nComponents = nCh;
104 }
105
106 //----------------------------------------------------------------------------------------------------------
107 // 1. Center: subtract per-channel mean
108 //----------------------------------------------------------------------------------------------------------
109 VectorXd vecMean = matData.rowwise().mean();
110 MatrixXd matCentered = matData.colwise() - vecMean;
111
112 //----------------------------------------------------------------------------------------------------------
113 // 2. Whiten
114 //----------------------------------------------------------------------------------------------------------
115 MatrixXd matWhitening, matDewhitening;
116 MatrixXd matWhite = whiten(matCentered, nComponents, matWhitening, matDewhitening);
117
118 //----------------------------------------------------------------------------------------------------------
119 // 3. FastICA — deflationary approach with logcosh (tanh) nonlinearity
120 //
121 // For each component i, iterate:
122 // g = tanh(W[i]^T * X_white) [1 x nSamples]
123 // g' = 1 - g^2 [1 x nSamples]
124 // w_new = (1/n) * X_white * g^T - mean(g') * w_i
125 // Gram-Schmidt orthogonalise against already-found rows
126 // Normalise
127 // Convergence: |w_new · w_old| ≥ 1 - tol
128 //----------------------------------------------------------------------------------------------------------
129 std::mt19937 rng(static_cast<unsigned>(randomSeed));
130 std::normal_distribution<double> dist(0.0, 1.0);
131
132 MatrixXd W_ica(nComponents, nComponents); // unmixing in whitened space
133 bool bConverged = true;
134
135 for (int i = 0; i < nComponents; ++i) {
136 // Random unit-vector initialisation
137 VectorXd w(nComponents);
138 for (int k = 0; k < nComponents; ++k) {
139 w(k) = dist(rng);
140 }
141 w.normalize();
142
143 bool compConverged = false;
144 for (int iter = 0; iter < maxIter; ++iter) {
145 // g = tanh(w^T * X_white) → row vector (1 x nSamples)
146 RowVectorXd u = w.transpose() * matWhite;
147 RowVectorXd g = u.array().tanh();
148 RowVectorXd gPrime = 1.0 - g.array().square(); // derivative of tanh
149
150 // Newton update
151 VectorXd wNew = (matWhite * g.transpose()) / nSamples
152 - gPrime.mean() * w;
153
154 // Deflation: orthogonalise against already-converged components
155 gramSchmidt(wNew, W_ica, i);
156
157 // Normalise
158 double norm = wNew.norm();
159 if (norm < 1e-12) {
160 // Degenerate — reinitialise randomly and restart this component
161 for (int k = 0; k < nComponents; ++k) {
162 w(k) = dist(rng);
163 }
164 w.normalize();
165 continue;
166 }
167 wNew /= norm;
168
169 // Convergence check: |w_new · w_old| should approach 1
170 double delta = std::abs(wNew.dot(w)) - 1.0;
171 w = wNew;
172
173 if (std::abs(delta) < tol) {
174 compConverged = true;
175 break;
176 }
177 }
178
179 if (!compConverged) {
180 bConverged = false;
181 qWarning() << "ICA::run: component" << i << "did not converge within" << maxIter << "iterations.";
182 }
183
184 W_ica.row(i) = w.transpose();
185 }
186
187 //----------------------------------------------------------------------------------------------------------
188 // 4. Compose full (sensor-space) unmixing and mixing matrices
189 //
190 // W_full = W_ica * W_whitening (n_comp x n_ch)
191 // A_full = W_dewhitening * W_ica^T (n_ch x n_comp) — exact inverse when n_comp == n_ch,
192 // pseudo-inverse otherwise
193 //----------------------------------------------------------------------------------------------------------
194 MatrixXd matUnmixing = W_ica * matWhitening; // n_comp x n_ch
195 MatrixXd matMixing = matDewhitening * W_ica.transpose(); // n_ch x n_comp
196
197 //----------------------------------------------------------------------------------------------------------
198 // 5. Compute source time series
199 //----------------------------------------------------------------------------------------------------------
200 MatrixXd matSources = matUnmixing * matCentered; // n_comp x n_samples
201
202 IcaResult result;
203 result.matMixing = std::move(matMixing);
204 result.matUnmixing = std::move(matUnmixing);
205 result.matSources = std::move(matSources);
206 result.vecMean = std::move(vecMean);
207 result.bConverged = bConverged;
208
209 return result;
210}
211
212//=============================================================================================================
213
214MatrixXd ICA::applyUnmixing(const MatrixXd& matData, const IcaResult& result)
215{
216 // Centre using the training mean, then project
217 MatrixXd matCentered = matData.colwise() - result.vecMean;
218 return result.matUnmixing * matCentered;
219}
220
221//=============================================================================================================
222
223MatrixXd ICA::excludeComponents(const MatrixXd& matData,
224 const IcaResult& result,
225 const QVector<int>& excludedComponents)
226{
227 if (excludedComponents.isEmpty()) {
228 return matData;
229 }
230
231 // Project to source space
232 MatrixXd matSources = applyUnmixing(matData, result);
233
234 // Zero the excluded components
235 for (int idx : excludedComponents) {
236 if (idx >= 0 && idx < static_cast<int>(matSources.rows())) {
237 matSources.row(idx).setZero();
238 }
239 }
240
241 // Back-project and add mean
242 MatrixXd matRecon = result.matMixing * matSources;
243 return matRecon.colwise() + result.vecMean;
244}
245
246//=============================================================================================================
247
248MatrixXd ICA::whiten(const MatrixXd& matCentered,
249 int nComponents,
250 MatrixXd& matWhitening,
251 MatrixXd& matDewhitening)
252{
253 const int nSamples = static_cast<int>(matCentered.cols());
254
255 // Sample covariance (unscaled)
256 MatrixXd cov = matCentered * matCentered.transpose() / static_cast<double>(nSamples);
257
258 // Eigendecomposition — SelfAdjointEigenSolver returns eigenvalues in ascending order
259 SelfAdjointEigenSolver<MatrixXd> eig(cov);
260
261 const int nCh = static_cast<int>(cov.rows());
262 // Take the nComponents largest eigenvalues/vectors (rightmost columns)
263 VectorXd eigenvalues = eig.eigenvalues().tail(nComponents).cwiseMax(1e-12);
264 MatrixXd eigenvectors = eig.eigenvectors().rightCols(nComponents); // n_ch x n_comp
265
266 // Whitening : W_w = D^{-1/2} * V^T (n_comp x n_ch)
267 // Dewhitening: W_d = V * D^{1/2} (n_ch x n_comp)
268 VectorXd sqrtEig = eigenvalues.cwiseSqrt();
269 VectorXd invSqrtEig = sqrtEig.cwiseInverse();
270
271 matWhitening = invSqrtEig.asDiagonal() * eigenvectors.transpose();
272 matDewhitening = eigenvectors * sqrtEig.asDiagonal();
273
274 return matWhitening * matCentered;
275}
Declaration of the ICA class providing FastICA-based independent component analysis for artifact remo...
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
Result of an ICA decomposition.
Definition ica.h:76
Eigen::MatrixXd matUnmixing
Definition ica.h:78
Eigen::MatrixXd matSources
Definition ica.h:79
Eigen::MatrixXd matMixing
Definition ica.h:77
Eigen::VectorXd vecMean
Definition ica.h:80
bool bConverged
Definition ica.h:81
static Eigen::MatrixXd excludeComponents(const Eigen::MatrixXd &matData, const IcaResult &result, const QVector< int > &excludedComponents)
Definition ica.cpp:223
static IcaResult run(const Eigen::MatrixXd &matData, int nComponents=-1, int maxIter=200, double tol=1e-4, int randomSeed=42)
Definition ica.cpp:93
static Eigen::MatrixXd applyUnmixing(const Eigen::MatrixXd &matData, const IcaResult &result)
Definition ica.cpp:214