44 #include "../network/networknode.h"
45 #include "../network/networkedge.h"
46 #include "../network/network.h"
55 #include <QtConcurrent>
61 #include <unsupported/Eigen/FFT>
67 using namespace CONNECTIVITYLIB;
68 using namespace Eigen;
69 using namespace UTILSLIB;
93 if(connectivitySettings.isEmpty()) {
94 qDebug() <<
"PhaseLockingValue::calculate - Input data is empty";
98 if(AbstractMetric::m_bStorageModeIsActive ==
false) {
99 connectivitySettings.clearIntermediateData();
104 #ifdef EIGEN_FFTW_DEFAULT
105 fftw_make_planner_thread_safe();
109 int iNRows = connectivitySettings.at(0).matData.rows();
110 RowVectorXf rowVert = RowVectorXf::Zero(3);
112 for(
int i = 0; i < iNRows; ++i) {
113 rowVert = RowVectorXf::Zero(3);
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);
125 int iSignalLength = connectivitySettings.at(0).matData.cols();
126 int iNfft = connectivitySettings.getFFTSize();
129 QPair<MatrixXd, VectorXd> tapers = Spectral::generateTapers(iSignalLength, connectivitySettings.getWindowType());
132 int iNFreqs = int(floor(iNfft / 2.0)) + 1;
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;
153 connectivitySettings.getIntermediateSumData().vecPairCsdSum,
154 connectivitySettings.getIntermediateSumData().vecPairCsdNormalizedSum,
167 QFuture<void> result = QtConcurrent::map(connectivitySettings.getTrialData(),
169 result.waitForFinished();
176 computePLV(connectivitySettings,
189 QVector<QPair<int,Eigen::MatrixXcd> >& vecPairCsdSum,
190 QVector<QPair<int,MatrixXcd> >& vecPairCsdNormalizedSum,
195 const QPair<MatrixXd, VectorXd>& tapers)
197 if(inputData.vecPairCsdNormalized.size() == iNRows) {
206 if(inputData.vecTapSpectra.isEmpty()) {
207 RowVectorXd vecInputFFT, rowData;
208 RowVectorXcd vecTmpFreq;
210 MatrixXcd matTapSpectrum(tapers.first.rows(), iNFreqs);
213 fft.SetFlag(fft.HalfSpectrum);
215 for (i = 0; i < iNRows; ++i) {
217 rowData.array() = inputData.matData.row(i).array() - inputData.matData.row(i).mean();
220 for(j = 0; j < tapers.first.rows(); j++) {
222 if (rowData.cols() < iNfft) {
223 vecInputFFT.setZero(iNfft);
224 vecInputFFT.block(0,0,1,rowData.cols()) = rowData.cwiseProduct(tapers.first.row(j));;
226 vecInputFFT = rowData.cwiseProduct(tapers.first.row(j));
230 fft.fwd(vecTmpFreq, vecInputFFT, iNfft);
231 matTapSpectrum.row(j) = vecTmpFreq * tapers.second(j);
234 inputData.vecTapSpectra.append(matTapSpectrum);
239 if(inputData.vecPairCsd.isEmpty()) {
240 MatrixXcd matCsd = MatrixXcd(iNRows, m_iNumberBinAmount);
242 bool bNfftEven =
false;
247 double denomCSD = sqrt(tapers.second.cwiseAbs2().sum()) * sqrt(tapers.second.cwiseAbs2().sum()) / 2.0;
249 for (i = 0; i < iNRows; ++i) {
250 for (j = i; j < iNRows; ++j) {
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;
255 if(m_iNumberBinStart == 0) {
256 matCsd.row(j)(0) /= 2.0;
259 if(bNfftEven && m_iNumberBinStart + m_iNumberBinAmount >= iNFreqs) {
260 matCsd.row(j).tail(1) /= 2.0;
264 inputData.vecPairCsd.append(QPair<int,MatrixXcd>(i,matCsd));
265 inputData.vecPairCsdNormalized.append(QPair<int,MatrixXcd>(i,matCsd.cwiseQuotient(matCsd.cwiseAbs())));
270 if(vecPairCsdSum.isEmpty()) {
271 vecPairCsdSum = inputData.vecPairCsd;
272 vecPairCsdNormalizedSum = inputData.vecPairCsdNormalized;
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;
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())));
289 if(vecPairCsdNormalizedSum.isEmpty()) {
290 vecPairCsdNormalizedSum = inputData.vecPairCsdNormalized;
292 for (
int j = 0; j < vecPairCsdNormalizedSum.size(); ++j) {
293 vecPairCsdNormalizedSum[j].second += inputData.vecPairCsdNormalized.at(j).second;
301 if(!m_bStorageModeIsActive) {
302 inputData.vecPairCsd.clear();
303 inputData.vecTapSpectra.clear();
304 inputData.vecPairCsdNormalized.clear();
316 QSharedPointer<NetworkEdge> pEdge;
319 for (
int i = 0; i < connectivitySettings.at(0).matData.rows(); ++i) {
320 matNom = connectivitySettings.getIntermediateSumData().vecPairCsdNormalizedSum.at(i).second.cwiseAbs() / connectivitySettings.size();
322 for(j = i; j < connectivitySettings.at(0).matData.rows(); ++j) {
323 matWeight = matNom.row(j).transpose();
325 pEdge = QSharedPointer<NetworkEdge>(
new NetworkEdge(i, j, matWeight));
327 finalNetwork.
getNodeAt(i)->append(pEdge);
328 finalNetwork.
getNodeAt(j)->append(pEdge);
329 finalNetwork.
append(pEdge);