MNE-CPP 0.1.9
A Framework for Electrophysiology
Loading...
Searching...
No Matches
phaselockingvalue.cpp
Go to the documentation of this file.
1//=============================================================================================================
39//=============================================================================================================
40// INCLUDES
41//=============================================================================================================
42
43#include "phaselockingvalue.h"
44#include "../network/networknode.h"
45#include "../network/networkedge.h"
46#include "../network/network.h"
47
48#include <utils/spectral.h>
49
50//=============================================================================================================
51// QT INCLUDES
52//=============================================================================================================
53
54#include <QDebug>
55#include <QtConcurrent>
56
57//=============================================================================================================
58// EIGEN INCLUDES
59//=============================================================================================================
60
61#include <unsupported/Eigen/FFT>
62
63//=============================================================================================================
64// USED NAMESPACES
65//=============================================================================================================
66
67using namespace CONNECTIVITYLIB;
68using namespace Eigen;
69using namespace UTILSLIB;
70
71//=============================================================================================================
72// DEFINE GLOBAL METHODS
73//=============================================================================================================
74
75//=============================================================================================================
76// DEFINE MEMBER METHODS
77//=============================================================================================================
78
82
83//*******************************************************************************************************
84
86{
87// QElapsedTimer timer;
88// qint64 iTime = 0;
89// timer.start();
90
91 Network finalNetwork("PLV");
92
93 if(connectivitySettings.isEmpty()) {
94 qDebug() << "PhaseLockingValue::calculate - Input data is empty";
95 return finalNetwork;
96 }
97
98 if(AbstractMetric::m_bStorageModeIsActive == false) {
99 connectivitySettings.clearIntermediateData();
100 }
101
102 finalNetwork.setSamplingFrequency(connectivitySettings.getSamplingFrequency());
103
104 #ifdef EIGEN_FFTW_DEFAULT
105 fftw_make_planner_thread_safe();
106 #endif
107
108 //Create nodes
109 int iNRows = connectivitySettings.at(0).matData.rows();
110 RowVectorXf rowVert = RowVectorXf::Zero(3);
111
112 for(int i = 0; i < iNRows; ++i) {
113 rowVert = RowVectorXf::Zero(3);
114
115 if(connectivitySettings.getNodePositions().rows() != 0 && i < connectivitySettings.getNodePositions().rows()) {
116 rowVert(0) = connectivitySettings.getNodePositions().row(i)(0);
117 rowVert(1) = connectivitySettings.getNodePositions().row(i)(1);
118 rowVert(2) = connectivitySettings.getNodePositions().row(i)(2);
119 }
120
121 finalNetwork.append(NetworkNode::SPtr(new NetworkNode(i, rowVert)));
122 }
123
124 // Check that iNfft >= signal length
125 int iSignalLength = connectivitySettings.at(0).matData.cols();
126 int iNfft = connectivitySettings.getFFTSize();
127
128 // Generate tapers
129 QPair<MatrixXd, VectorXd> tapers = Spectral::generateTapers(iSignalLength, connectivitySettings.getWindowType());
130
131 // Initialize
132 int iNFreqs = int(floor(iNfft / 2.0)) + 1;
133
134 // Check if start and bin amount need to be reset to full spectrum
135 if(m_iNumberBinStart == -1 ||
136 m_iNumberBinAmount == -1 ||
137 m_iNumberBinStart > iNFreqs ||
138 m_iNumberBinAmount > iNFreqs ||
139 m_iNumberBinAmount + m_iNumberBinStart > iNFreqs) {
140 qDebug() << "PhaseLockingValue::calculate - Resetting to full spectrum";
141 AbstractMetric::m_iNumberBinStart = 0;
142 AbstractMetric::m_iNumberBinAmount = iNFreqs;
143 }
144
145 // Pass information about the FFT length. Use iNFreqs because we only use the half spectrum
146 finalNetwork.setFFTSize(iNFreqs);
147 finalNetwork.setUsedFreqBins(AbstractMetric::m_iNumberBinAmount);
148
149 QMutex mutex;
150
151 std::function<void(ConnectivitySettings::IntermediateTrialData&)> computeLambda = [&](ConnectivitySettings::IntermediateTrialData& inputData) {
152 compute(inputData,
153 connectivitySettings.getIntermediateSumData().vecPairCsdSum,
154 connectivitySettings.getIntermediateSumData().vecPairCsdNormalizedSum,
155 mutex,
156 iNRows,
157 iNFreqs,
158 iNfft,
159 tapers);
160 };
161
162// iTime = timer.elapsed();
163// qWarning() << "Preparation" << iTime;
164// timer.restart();
165
166 // Compute PLV in parallel for all trials
167 QFuture<void> result = QtConcurrent::map(connectivitySettings.getTrialData(),
168 computeLambda);
169 result.waitForFinished();
170
171// iTime = timer.elapsed();
172// qWarning() << "ComputeSpectraPSDCSD" << iTime;
173// timer.restart();
174
175 // Compute PLV
176 computePLV(connectivitySettings,
177 finalNetwork);
178
179// iTime = timer.elapsed();
180// qWarning() << "Compute" << iTime;
181// timer.restart();
182
183 return finalNetwork;
184}
185
186//=============================================================================================================
187
189 QVector<QPair<int,Eigen::MatrixXcd> >& vecPairCsdSum,
190 QVector<QPair<int,MatrixXcd> >& vecPairCsdNormalizedSum,
191 QMutex& mutex,
192 int iNRows,
193 int iNFreqs,
194 int iNfft,
195 const QPair<MatrixXd, VectorXd>& tapers)
196{
197 if(inputData.vecPairCsdNormalized.size() == iNRows) {
198 //qDebug() << "PhaseLockingValue::compute - vecPairCsdNormalized was already computed for this trial.";
199 return;
200 }
201
202 int i,j;
203
204 // Calculate tapered spectra if not available already
205 // This code was copied and changed modified Utils/Spectra since we do not want to call the function due to time loss.
206 if(inputData.vecTapSpectra.isEmpty()) {
207 RowVectorXd vecInputFFT, rowData;
208 RowVectorXcd vecTmpFreq;
209
210 MatrixXcd matTapSpectrum(tapers.first.rows(), iNFreqs);
211
212 FFT<double> fft;
213 fft.SetFlag(fft.HalfSpectrum);
214
215 for (i = 0; i < iNRows; ++i) {
216 // Substract mean
217 rowData.array() = inputData.matData.row(i).array() - inputData.matData.row(i).mean();
218
219 // Calculate tapered spectra if not available already
220 for(j = 0; j < tapers.first.rows(); j++) {
221 // Zero padd if necessary. The zero padding in Eigen's FFT is only working for column vectors.
222 if (rowData.cols() < iNfft) {
223 vecInputFFT.setZero(iNfft);
224 vecInputFFT.block(0,0,1,rowData.cols()) = rowData.cwiseProduct(tapers.first.row(j));;
225 } else {
226 vecInputFFT = rowData.cwiseProduct(tapers.first.row(j));
227 }
228
229 // FFT for freq domain returning the half spectrum and multiply taper weights
230 fft.fwd(vecTmpFreq, vecInputFFT, iNfft);
231 matTapSpectrum.row(j) = vecTmpFreq * tapers.second(j);
232 }
233
234 inputData.vecTapSpectra.append(matTapSpectrum);
235 }
236 }
237
238 // Compute CSD
239 if(inputData.vecPairCsd.isEmpty()) {
240 MatrixXcd matCsd = MatrixXcd(iNRows, m_iNumberBinAmount);
241
242 bool bNfftEven = false;
243 if (iNfft % 2 == 0){
244 bNfftEven = true;
245 }
246
247 double denomCSD = sqrt(tapers.second.cwiseAbs2().sum()) * sqrt(tapers.second.cwiseAbs2().sum()) / 2.0;
248
249 for (i = 0; i < iNRows; ++i) {
250 for (j = i; j < iNRows; ++j) {
251 // Compute CSD (average over tapers if necessary)
252 matCsd.row(j) = inputData.vecTapSpectra.at(i).block(0,m_iNumberBinStart,inputData.vecTapSpectra.at(i).rows(),m_iNumberBinAmount).cwiseProduct(inputData.vecTapSpectra.at(j).block(0,m_iNumberBinStart,inputData.vecTapSpectra.at(j).rows(),m_iNumberBinAmount).conjugate()).colwise().sum() / denomCSD;
253
254 // Divide first and last element by 2 due to half spectrum
255 if(m_iNumberBinStart == 0) {
256 matCsd.row(j)(0) /= 2.0;
257 }
258
259 if(bNfftEven && m_iNumberBinStart + m_iNumberBinAmount >= iNFreqs) {
260 matCsd.row(j).tail(1) /= 2.0;
261 }
262 }
263
264 inputData.vecPairCsd.append(QPair<int,MatrixXcd>(i,matCsd));
265 inputData.vecPairCsdNormalized.append(QPair<int,MatrixXcd>(i,matCsd.cwiseQuotient(matCsd.cwiseAbs())));
266 }
267
268 mutex.lock();
269
270 if(vecPairCsdSum.isEmpty()) {
271 vecPairCsdSum = inputData.vecPairCsd;
272 vecPairCsdNormalizedSum = inputData.vecPairCsdNormalized;
273 } else {
274 for (int j = 0; j < vecPairCsdSum.size(); ++j) {
275 vecPairCsdSum[j].second += inputData.vecPairCsd.at(j).second;
276 vecPairCsdNormalizedSum[j].second += inputData.vecPairCsdNormalized.at(j).second;
277 }
278 }
279
280 mutex.unlock();
281 } else {
282 if(inputData.vecPairCsdNormalized.isEmpty()) {
283 for (i = 0; i < iNRows; ++i) {
284 inputData.vecPairCsdNormalized.append(QPair<int,MatrixXcd>(i,inputData.vecPairCsd.at(i).second.cwiseQuotient(inputData.vecPairCsd.at(i).second.cwiseAbs())));
285 }
286
287 mutex.lock();
288
289 if(vecPairCsdNormalizedSum.isEmpty()) {
290 vecPairCsdNormalizedSum = inputData.vecPairCsdNormalized;
291 } else {
292 for (int j = 0; j < vecPairCsdNormalizedSum.size(); ++j) {
293 vecPairCsdNormalizedSum[j].second += inputData.vecPairCsdNormalized.at(j).second;
294 }
295 }
296
297 mutex.unlock();
298 }
299 }
300
301 if(!m_bStorageModeIsActive) {
302 inputData.vecPairCsd.clear();
303 inputData.vecTapSpectra.clear();
304 inputData.vecPairCsdNormalized.clear();
305 }
306}
307
308//=============================================================================================================
309
311 Network& finalNetwork)
312{
313 // Compute final PLV and create Network
314 MatrixXd matNom;
315 MatrixXd matWeight;
316 QSharedPointer<NetworkEdge> pEdge;
317 int j;
318
319 for (int i = 0; i < connectivitySettings.at(0).matData.rows(); ++i) {
320 matNom = connectivitySettings.getIntermediateSumData().vecPairCsdNormalizedSum.at(i).second.cwiseAbs() / connectivitySettings.size();
321
322 for(j = i; j < connectivitySettings.at(0).matData.rows(); ++j) {
323 matWeight = matNom.row(j).transpose();
324
325 pEdge = QSharedPointer<NetworkEdge>(new NetworkEdge(i, j, matWeight));
326
327 finalNetwork.getNodeAt(i)->append(pEdge);
328 finalNetwork.getNodeAt(j)->append(pEdge);
329 finalNetwork.append(pEdge);
330 }
331 }
332}
PhaseLockingValue class declaration.
Declaration of Spectral class.
This class is a container for connectivity settings.
static void computePLV(ConnectivitySettings &connectivitySettings, Network &finalNetwork)
static void compute(ConnectivitySettings::IntermediateTrialData &inputData, QVector< QPair< int, Eigen::MatrixXcd > > &vecPairCsdSum, QVector< QPair< int, Eigen::MatrixXcd > > &vecPairCsdNormalizedSum, QMutex &mutex, int iNRows, int iNFreqs, int iNfft, const QPair< Eigen::MatrixXd, Eigen::VectorXd > &tapers)
static Network calculate(ConnectivitySettings &connectivitySettings)
This class holds information about a network, can compute a distance table and provide network metric...
Definition network.h:89
void setUsedFreqBins(int iNumberFreqBins)
Definition network.cpp:506
void append(QSharedPointer< NetworkEdge > newEdge)
void setFFTSize(int iFFTSize)
Definition network.cpp:513
void setSamplingFrequency(float fSFreq)
Definition network.cpp:492
QSharedPointer< NetworkNode > getNodeAt(int i)
Definition network.cpp:163
This class holds an object to describe the edge of a network.
Definition networkedge.h:80
This class holds an object to describe the node of a network.
Definition networknode.h:82
QSharedPointer< NetworkNode > SPtr
Definition networknode.h:85
static QPair< Eigen::MatrixXd, Eigen::VectorXd > generateTapers(int iSignalLength, const QString &sWindowType="hanning")
Definition spectral.cpp:292