v2.0.0
Loading...
Searching...
No Matches
mvar_model.cpp
Go to the documentation of this file.
1//=============================================================================================================
34
35//=============================================================================================================
36// INCLUDES
37//=============================================================================================================
38
39#include "mvar_model.h"
40
41//=============================================================================================================
42// QT INCLUDES
43//=============================================================================================================
44
45#include <QDebug>
46#include <QtMath>
47
48//=============================================================================================================
49// EIGEN INCLUDES
50//=============================================================================================================
51
52#include <Eigen/Dense>
53
54//=============================================================================================================
55// USED NAMESPACES
56//=============================================================================================================
57
58using namespace CONNLIB;
59using namespace Eigen;
60
61//=============================================================================================================
62// DEFINE GLOBAL METHODS
63//=============================================================================================================
64
65//=============================================================================================================
66// DEFINE MEMBER METHODS
67//=============================================================================================================
68
69void MvarModel::fit(const MatrixXd& data, int p)
70{
71 m_nChannels = static_cast<int>(data.rows());
72
73 if(p <= 0) {
74 p = selectOrderBIC(data);
75 }
76
77 fitLevinsonDurbin(data, p);
78}
79
80//=============================================================================================================
81
82QVector<MatrixXd> MvarModel::coefficients() const
83{
84 return m_coeffs;
85}
86
87//=============================================================================================================
88
89MatrixXd MvarModel::noiseCov() const
90{
91 return m_noiseCov;
92}
93
94//=============================================================================================================
95
97{
98 return m_order;
99}
100
101//=============================================================================================================
102
103QVector<MatrixXcd> MvarModel::transferFunction(const VectorXd& freqs) const
104{
105 QVector<MatrixXcd> vecH;
106 vecH.reserve(static_cast<int>(freqs.size()));
107
108 const int nCh = m_nChannels;
109 const MatrixXcd matI = MatrixXcd::Identity(nCh, nCh);
110 const std::complex<double> j(0.0, 1.0);
111
112 for(int fi = 0; fi < freqs.size(); ++fi) {
113 MatrixXcd matA = matI;
114
115 for(int k = 0; k < m_order; ++k) {
116 const double phase = -2.0 * M_PI * freqs(fi) * (k + 1);
117 const std::complex<double> expVal = std::exp(j * phase);
118 matA -= m_coeffs[k].cast<std::complex<double>>() * expVal;
119 }
120
121 vecH.append(matA.inverse());
122 }
123
124 return vecH;
125}
126
127//=============================================================================================================
128
129QVector<MatrixXcd> MvarModel::spectralMatrix(const VectorXd& freqs) const
130{
131 QVector<MatrixXcd> vecH = transferFunction(freqs);
132 QVector<MatrixXcd> vecS;
133 vecS.reserve(vecH.size());
134
135 const MatrixXcd matSigma = m_noiseCov.cast<std::complex<double>>();
136
137 for(int fi = 0; fi < vecH.size(); ++fi) {
138 vecS.append(vecH[fi] * matSigma * vecH[fi].adjoint());
139 }
140
141 return vecS;
142}
143
144//=============================================================================================================
145
146void MvarModel::fitLevinsonDurbin(const MatrixXd& data, int p)
147{
148 m_order = p;
149 m_nChannels = static_cast<int>(data.rows());
150
151 const int nCh = m_nChannels;
152 const int nSamples = static_cast<int>(data.cols());
153 const int nObs = nSamples - p;
154
155 if(nObs <= 0) {
156 qWarning() << "MvarModel::fitLevinsonDurbin - Not enough samples for model order" << p;
157 m_coeffs.clear();
158 m_noiseCov = MatrixXd::Identity(nCh, nCh);
159 return;
160 }
161
162 // Subtract mean from each channel
163 MatrixXd dataCentered = data;
164 for(int i = 0; i < nCh; ++i) {
165 dataCentered.row(i).array() -= dataCentered.row(i).mean();
166 }
167
168 // Build regression matrices for OLS solution of Yule-Walker equations
169 // Y = data(:, p:end) -> nCh x nObs
170 // Z = [data(:, p-1:end-1); -> nCh*p x nObs
171 // data(:, p-2:end-2);
172 // ...
173 // data(:, 0:end-p)]
174 MatrixXd matY = dataCentered.rightCols(nObs);
175 MatrixXd matZ(nCh * p, nObs);
176
177 for(int k = 0; k < p; ++k) {
178 matZ.middleRows(k * nCh, nCh) = dataCentered.middleCols(p - 1 - k, nObs);
179 }
180
181 // Solve Y = A * Z via least squares: A = Y * Z^T * (Z * Z^T)^{-1}
182 MatrixXd matZZT = matZ * matZ.transpose();
183 MatrixXd matYZT = matY * matZ.transpose();
184 MatrixXd matA = matYZT * matZZT.ldlt().solve(MatrixXd::Identity(nCh * p, nCh * p));
185
186 // Extract coefficient matrices A_1..A_p
187 m_coeffs.clear();
188 m_coeffs.reserve(p);
189 for(int k = 0; k < p; ++k) {
190 m_coeffs.append(matA.middleCols(k * nCh, nCh));
191 }
192
193 // Compute noise covariance from residuals
194 MatrixXd matE = matY - matA * matZ;
195 m_noiseCov = (matE * matE.transpose()) / static_cast<double>(nObs);
196}
197
198//=============================================================================================================
199
200int MvarModel::selectOrderBIC(const MatrixXd& data, int maxOrder) const
201{
202 const int nCh = static_cast<int>(data.rows());
203 const int nSamples = static_cast<int>(data.cols());
204
205 // Limit max order to avoid underdetermined systems
206 maxOrder = qMin(maxOrder, nSamples / (nCh + 1));
207 if(maxOrder < 1) {
208 maxOrder = 1;
209 }
210
211 // Subtract mean
212 MatrixXd dataCentered = data;
213 for(int i = 0; i < nCh; ++i) {
214 dataCentered.row(i).array() -= dataCentered.row(i).mean();
215 }
216
217 double bestBIC = std::numeric_limits<double>::max();
218 int bestOrder = 1;
219
220 for(int p = 1; p <= maxOrder; ++p) {
221 const int nObs = nSamples - p;
222 if(nObs <= nCh * p) {
223 break;
224 }
225
226 // Build regression matrices
227 MatrixXd matY = dataCentered.rightCols(nObs);
228 MatrixXd matZ(nCh * p, nObs);
229
230 for(int k = 0; k < p; ++k) {
231 matZ.middleRows(k * nCh, nCh) = dataCentered.middleCols(p - 1 - k, nObs);
232 }
233
234 // Solve via OLS
235 MatrixXd matZZT = matZ * matZ.transpose();
236 MatrixXd matYZT = matY * matZ.transpose();
237 MatrixXd matA = matYZT * matZZT.ldlt().solve(MatrixXd::Identity(nCh * p, nCh * p));
238
239 // Residual covariance
240 MatrixXd matE = matY - matA * matZ;
241 MatrixXd matSigma = (matE * matE.transpose()) / static_cast<double>(nObs);
242
243 // BIC = n * ln(det(Sigma)) + k * ln(n), where k = p * nCh^2
244 double detSigma = matSigma.determinant();
245 if(detSigma <= 0.0) {
246 continue;
247 }
248
249 const int nParams = p * nCh * nCh;
250 double bic = nObs * std::log(detSigma) + nParams * std::log(static_cast<double>(nObs));
251
252 if(bic < bestBIC) {
253 bestBIC = bic;
254 bestOrder = p;
255 }
256 }
257
258 return bestOrder;
259}
#define M_PI
MvarModel class declaration.
Functional connectivity metrics (coherence, PLV, cross-correlation, etc.).
Eigen::MatrixXd noiseCov() const
QVector< Eigen::MatrixXcd > spectralMatrix(const Eigen::VectorXd &freqs) const
QVector< Eigen::MatrixXcd > transferFunction(const Eigen::VectorXd &freqs) const
int order() const
void fit(const Eigen::MatrixXd &data, int p=0)
QVector< Eigen::MatrixXd > coefficients() const