MNE-CPP  0.1.9
A Framework for Electrophysiology
coherency.cpp
Go to the documentation of this file.
1 //=============================================================================================================
39 //=============================================================================================================
40 // INCLUDES
41 //=============================================================================================================
42 
43 #include "coherency.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 
67 using namespace CONNECTIVITYLIB;
68 using namespace Eigen;
69 using namespace UTILSLIB;
70 
71 //=============================================================================================================
72 // DEFINE GLOBAL METHODS
73 //=============================================================================================================
74 
75 //=============================================================================================================
76 // DEFINE MEMBER METHODS
77 //=============================================================================================================
78 
80 {
81 }
82 
83 //=============================================================================================================
84 
85 void Coherency::calculateAbs(Network& finalNetwork,
86  ConnectivitySettings &connectivitySettings)
87 {
88 // QElapsedTimer timer;
89 // qint64 iTime = 0;
90 // timer.start();
91 
92  if(connectivitySettings.isEmpty()) {
93  qDebug() << "Coherency::calculateReal - Input data is empty";
94  return;
95  }
96 
97  #ifdef EIGEN_FFTW_DEFAULT
98  fftw_make_planner_thread_safe();
99  #endif
100 
101  int iSignalLength = connectivitySettings.at(0).matData.cols();
102  int iNfft = connectivitySettings.getFFTSize();
103 
104  // Generate tapers
105  QPair<MatrixXd, VectorXd> tapers = Spectral::generateTapers(iSignalLength, connectivitySettings.getWindowType());
106 
107  // Initialize vecPsdAvg and vecCsdAvg
108  int iNRows = connectivitySettings.at(0).matData.rows();
109  int iNFreqs = int(floor(iNfft / 2.0)) + 1;
110 
111  // Compute PSD/CSD for each trial
112  QMutex mutex;
113 
114  std::function<void(ConnectivitySettings::IntermediateTrialData&)> computeLambda = [&](ConnectivitySettings::IntermediateTrialData& inputData) {
115  compute(inputData,
116  connectivitySettings.getIntermediateSumData().matPsdSum,
117  connectivitySettings.getIntermediateSumData().vecPairCsdSum,
118  mutex,
119  iNRows,
120  iNFreqs,
121  iNfft,
122  tapers);
123  };
124 
125 // iTime = timer.elapsed();
126 // qWarning() << "Preparation" << iTime;
127 // timer.restart();
128 
129  QFuture<void> result = QtConcurrent::map(connectivitySettings.getTrialData(),
130  computeLambda);
131  result.waitForFinished();
132 
133 // iTime = timer.elapsed();
134 // qWarning() << "ComputeSpectraPSDCSD" << iTime;
135 // timer.restart();
136 
137  // Compute CSD/sqrt(PSD_X * PSD_Y)
138  std::function<void(QPair<int,MatrixXcd>&)> computePSDCSDLambda = [&](QPair<int,MatrixXcd>& pairInput) {
139  computePSDCSDAbs(mutex,
140  finalNetwork,
141  pairInput,
142  connectivitySettings.getIntermediateSumData().matPsdSum);
143  };
144 
145  QFuture<void> resultCSDPSD = QtConcurrent::map(connectivitySettings.getIntermediateSumData().vecPairCsdSum,
146  computePSDCSDLambda);
147  resultCSDPSD.waitForFinished();
148 
149 // iTime = timer.elapsed();
150 // qWarning() << "Compute" << iTime;
151 // timer.restart();
152 }
153 
154 //=============================================================================================================
155 
157  ConnectivitySettings &connectivitySettings)
158 {
159 // QElapsedTimer timer;
160 // qint64 iTime = 0;
161 // timer.start();
162 
163  if(connectivitySettings.isEmpty()) {
164  qDebug() << "Coherency::calculateImag - Input data is empty";
165  return;
166  }
167 
168  #ifdef EIGEN_FFTW_DEFAULT
169  fftw_make_planner_thread_safe();
170  #endif
171 
172  int iSignalLength = connectivitySettings.at(0).matData.cols();
173  int iNfft = connectivitySettings.getFFTSize();
174 
175  // Generate tapers
176  QPair<MatrixXd, VectorXd> tapers = Spectral::generateTapers(iSignalLength, connectivitySettings.getWindowType());
177 
178  // Initialize vecPsdAvg and vecCsdAvg
179  int iNRows = connectivitySettings.at(0).matData.rows();
180  int iNFreqs = int(floor(iNfft / 2.0)) + 1;
181 
182  // Compute PSD/CSD for each trial
183  QMutex mutex;
184 
185  std::function<void(ConnectivitySettings::IntermediateTrialData&)> computeLambda = [&](ConnectivitySettings::IntermediateTrialData& inputData) {
186  compute(inputData,
187  connectivitySettings.getIntermediateSumData().matPsdSum,
188  connectivitySettings.getIntermediateSumData().vecPairCsdSum,
189  mutex,
190  iNRows,
191  iNFreqs,
192  iNfft,
193  tapers);
194  };
195 
196 // iTime = timer.elapsed();
197 // qWarning() << "Preparation" << iTime;
198 // timer.restart();
199 
200  QFuture<void> result = QtConcurrent::map(connectivitySettings.getTrialData(),
201  computeLambda);
202  result.waitForFinished();
203 
204 // iTime = timer.elapsed();
205 // qWarning() << "ComputeSpectraPSDCSD" << iTime;
206 // timer.restart();
207 
208  // Compute CSD/sqrt(PSD_X * PSD_Y)
209  std::function<void(QPair<int,MatrixXcd>&)> computePSDCSDLambda = [&](QPair<int,MatrixXcd>& pairInput) {
210  computePSDCSDImag(mutex,
211  finalNetwork,
212  pairInput,
213  connectivitySettings.getIntermediateSumData().matPsdSum);
214  };
215 
216  QFuture<void> resultCSDPSD = QtConcurrent::map(connectivitySettings.getIntermediateSumData().vecPairCsdSum,
217  computePSDCSDLambda);
218  resultCSDPSD.waitForFinished();
219 
220 // iTime = timer.elapsed();
221 // qWarning() << "Compute" << iTime;
222 // timer.restart();
223 }
224 
225 //=============================================================================================================
226 
227 void Coherency::compute(ConnectivitySettings::IntermediateTrialData& inputData,
228  MatrixXd& matPsdSum,
229  QVector<QPair<int,MatrixXcd> >& vecPairCsdSum,
230  QMutex& mutex,
231  int iNRows,
232  int iNFreqs,
233  int iNfft,
234  const QPair<MatrixXd, VectorXd>& tapers)
235 {
236 // QElapsedTimer timer;
237 // qint64 iTime = 0;
238 // timer.start();
239 
240  if(inputData.vecPairCsd.size() == iNRows) {
241  //qDebug() << "Coherency::compute - vecPairCsd were already computed for this trial.";
242  return;
243  }
244 
245  //qDebug() << "Coherency::compute - vecPairCsdSum and matPsdSum are computed for this trial.";
246 
247  // Substract mean, compute tapered spectra and PSD
248  // This code was copied and changed modified Utils/Spectra since we do not want to call the function due to time loss.
249  bool bNfftEven = false;
250  if (iNfft % 2 == 0){
251  bNfftEven = true;
252  }
253 
254  FFT<double> fft;
255  fft.SetFlag(fft.HalfSpectrum);
256 
257  double denomPSD = tapers.second.cwiseAbs2().sum() / 2.0;
258 
259  RowVectorXd vecInputFFT, rowData;
260  RowVectorXcd vecTmpFreq;
261 
262  MatrixXcd matTapSpectrum(tapers.first.rows(), iNFreqs);
263 
264  int i,j;
265 
266  inputData.matPsd = MatrixXd(iNRows, m_iNumberBinAmount);
267 
268  for (i = 0; i < iNRows; ++i) {
269  // Substract mean
270  rowData.array() = inputData.matData.row(i).array() - inputData.matData.row(i).mean();
271 
272  // Calculate tapered spectra if not available already
273  if(inputData.vecTapSpectra.size() != iNRows) {
274  for(j = 0; j < tapers.first.rows(); j++) {
275  // Zero padd if necessary. The zero padding in Eigen's FFT is only working for column vectors.
276  if (rowData.cols() < iNfft) {
277  vecInputFFT.setZero(iNfft);
278  vecInputFFT.block(0,0,1,rowData.cols()) = rowData.cwiseProduct(tapers.first.row(j));;
279  } else {
280  vecInputFFT = rowData.cwiseProduct(tapers.first.row(j));
281  }
282 
283  // FFT for freq domain returning the half spectrum and multiply taper weights
284  fft.fwd(vecTmpFreq, vecInputFFT, iNfft);
285  matTapSpectrum.row(j) = vecTmpFreq * tapers.second(j);
286  }
287 
288  inputData.vecTapSpectra.append(matTapSpectrum);
289  }
290 
291  // Compute PSD (average over tapers if necessary).
292  inputData.matPsd.row(i) = inputData.vecTapSpectra.at(i).block(0,m_iNumberBinStart,inputData.vecTapSpectra.at(i).rows(),m_iNumberBinAmount).cwiseAbs2().colwise().sum() / denomPSD;
293 
294  // Divide first and last element by 2 due to half spectrum
295  if(m_iNumberBinStart == 0) {
296  inputData.matPsd.row(i)(0) /= 2.0;
297  }
298 
299  if(bNfftEven && m_iNumberBinStart + m_iNumberBinAmount >= iNFreqs) {
300  inputData.matPsd.row(i).tail(1) /= 2.0;
301  }
302  }
303 
304  mutex.lock();
305 
306  if(matPsdSum.rows() == 0 || matPsdSum.cols() == 0) {
307  matPsdSum = inputData.matPsd;
308  } else {
309  matPsdSum += inputData.matPsd;
310  }
311 
312  mutex.unlock();
313 
314 // iTime = timer.elapsed();
315 // qWarning() << QThread::currentThreadId() << "Coherency::compute timer - compute - Tapered spectra and PSD (summing):" << iTime;
316 // timer.restart();
317 
318  // Compute CSD
319  if(inputData.vecPairCsd.size() != iNRows) {
320  inputData.vecPairCsd.clear();
321 
322  //MatrixXcd matCsd = MatrixXcd(iNRows, iNFreqs);
323  MatrixXcd matCsd = MatrixXcd(iNRows, m_iNumberBinAmount);
324 
325  double denomCSD = sqrt(tapers.second.cwiseAbs2().sum()) * sqrt(tapers.second.cwiseAbs2().sum()) / 2.0;
326 
327  for (i = 0; i < iNRows; ++i) {
328  for (j = i; j < iNRows; ++j) {
329  // Compute CSD (average over tapers if necessary)
330  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;
331 
332  // Divide first and last element by 2 due to half spectrum
333  if(m_iNumberBinStart == 0) {
334  matCsd.row(j)(0) /= 2.0;
335  }
336 
337  if(bNfftEven && m_iNumberBinStart + m_iNumberBinAmount >= iNFreqs) {
338  matCsd.row(j).tail(1) /= 2.0;
339  }
340  }
341 
342  inputData.vecPairCsd.append(QPair<int,MatrixXcd>(i,matCsd));
343  }
344 
345  mutex.lock();
346 
347  if(vecPairCsdSum.isEmpty()) {
348  vecPairCsdSum = inputData.vecPairCsd;
349  } else {
350  for (j = 0; j < vecPairCsdSum.size(); ++j) {
351  vecPairCsdSum[j].second += inputData.vecPairCsd.at(j).second;
352  }
353  }
354 
355  mutex.unlock();
356  }
357 
358 // iTime = timer.elapsed();
359 // qWarning() << QThread::currentThreadId() << "Coherency::compute timer - compute - CSD summing:" << iTime;
360 // timer.restart();
361 
362  //Do not store data to save memory
363  if(!m_bStorageModeIsActive) {
364  inputData.vecPairCsd.clear();
365  inputData.vecTapSpectra.clear();
366  }
367 
368 // iTime = timer.elapsed();
369 // qWarning() << QThread::currentThreadId() << "Coherency::compute timer - compute - Deleting data:" << iTime;
370 // timer.restart();
371 }
372 
373 //=============================================================================================================
374 
375 void Coherency::computePSDCSDAbs(QMutex& mutex,
376  Network& finalNetwork,
377  const QPair<int,MatrixXcd>& pairInput,
378  const MatrixXd& matPsdSum)
379 {
380  MatrixXd matPSDtmp(matPsdSum.rows(), matPsdSum.cols());
381  RowVectorXd rowPsdSum = matPsdSum.row(pairInput.first);
382 
383  for(int j = 0; j < matPSDtmp.rows(); ++j) {
384  matPSDtmp.row(j) = rowPsdSum.cwiseProduct(matPsdSum.row(j));
385  }
386 
387  // Average. Note that the number of trials cancel each other out.
388  MatrixXcd matCohy = pairInput.second.cwiseQuotient(matPSDtmp.cwiseSqrt());
389 
390  QSharedPointer<NetworkEdge> pEdge;
391  MatrixXd matWeight;
392  int j;
393  int i = pairInput.first;
394 
395  for(j = i; j < matCohy.rows(); ++j) {
396  matWeight = matCohy.row(j).cwiseAbs().transpose();
397  pEdge = QSharedPointer<NetworkEdge>(new NetworkEdge(i, j, matWeight));
398 
399  mutex.lock();
400  finalNetwork.getNodeAt(i)->append(pEdge);
401  finalNetwork.getNodeAt(j)->append(pEdge);
402  finalNetwork.append(pEdge);
403  mutex.unlock();
404  }
405 }
406 
407 //=============================================================================================================
408 
409 void Coherency::computePSDCSDImag(QMutex& mutex,
410  Network& finalNetwork,
411  const QPair<int,MatrixXcd>& pairInput,
412  const MatrixXd& matPsdSum)
413 {
414  MatrixXd matPSDtmp(matPsdSum.rows(), matPsdSum.cols());
415  RowVectorXd rowPsdSum = matPsdSum.row(pairInput.first);
416 
417  for(int j = 0; j < matPSDtmp.rows(); ++j) {
418  matPSDtmp.row(j) = rowPsdSum.cwiseProduct(matPsdSum.row(j));
419  }
420 
421  MatrixXcd matCohy = pairInput.second.cwiseQuotient(matPSDtmp.cwiseSqrt());
422 
423  QSharedPointer<NetworkEdge> pEdge;
424  MatrixXd matWeight;
425  int j;
426  int i = pairInput.first;
427 
428  for(j = i; j < matCohy.rows(); ++j) {
429  matWeight = matCohy.row(j).imag().transpose();
430  pEdge = QSharedPointer<NetworkEdge>(new NetworkEdge(i, j, matWeight));
431 
432  mutex.lock();
433  finalNetwork.getNodeAt(i)->append(pEdge);
434  finalNetwork.getNodeAt(j)->append(pEdge);
435  finalNetwork.append(pEdge);
436  mutex.unlock();
437  }
438 }
CONNECTIVITYLIB::Coherency::calculateAbs
static void calculateAbs(Network &finalNetwork, ConnectivitySettings &connectivitySettings)
Definition: coherency.cpp:85
CONNECTIVITYLIB::ConnectivitySettings::IntermediateTrialData
Definition: connectivitysettings.h:98
spectral.h
Declaration of Spectral class.
CONNECTIVITYLIB::Network
This class holds information about a network, can compute a distance table and provide network metric...
Definition: network.h:88
CONNECTIVITYLIB::Network::append
void append(QSharedPointer< NetworkEdge > newEdge)
CONNECTIVITYLIB::Network::getNodeAt
QSharedPointer< NetworkNode > getNodeAt(int i)
Definition: network.cpp:163
coherency.h
Coherency class declaration.
CONNECTIVITYLIB::Coherency::calculateImag
static void calculateImag(Network &finalNetwork, ConnectivitySettings &connectivitySettings)
Definition: coherency.cpp:156
CONNECTIVITYLIB::NetworkEdge
This class holds an object to describe the edge of a network.
Definition: networkedge.h:79
CONNECTIVITYLIB::ConnectivitySettings
This class is a container for connectivity settings.
Definition: connectivitysettings.h:91
CONNECTIVITYLIB::Coherency::Coherency
Coherency()
Definition: coherency.cpp:79