43#include <QtConcurrent>
44#include <QRandomGenerator>
63 const QVector<MatrixXd>& dataA,
64 const QVector<MatrixXd>& dataB,
65 const SparseMatrix<int>& adjacency,
73 const int nA = dataA.size();
76 MatrixXd tObs = computeTMap(dataA, dataB);
80 double threshold = inverseTCdf(clusterAlpha, df);
83 auto [clusterIds, clusterStats] = findClusters(tObs, threshold, adjacency, tail);
86 QVector<MatrixXd> allData;
87 allData.reserve(nA + dataB.size());
88 for (
const auto& m : dataA) allData.append(m);
89 for (
const auto& m : dataB) allData.append(m);
92 QVector<double> nullDist(nPermutations);
95 QVector<int> permIndices(nPermutations);
96 std::iota(permIndices.begin(), permIndices.end(), 0);
98 std::function<double(
int)> permuteFunc = [&](int) ->
double {
99 return permuteOnce(allData, nA, adjacency, threshold, tail);
102 QFuture<double> future = QtConcurrent::mapped(permIndices, permuteFunc);
103 future.waitForFinished();
105 for (
int i = 0; i < nPermutations; ++i) {
106 nullDist[i] = future.resultAt(i);
110 std::sort(nullDist.begin(), nullDist.end());
113 QVector<double> clusterPvals(clusterStats.size());
114 for (
int c = 0; c < clusterStats.size(); ++c) {
115 double obsStat = std::fabs(clusterStats[c]);
117 for (
int p = 0; p < nPermutations; ++p) {
118 if (nullDist[p] >= obsStat) {
122 clusterPvals[c] =
static_cast<double>(count) /
static_cast<double>(nPermutations);
136MatrixXd StatsCluster::computeTMap(
const QVector<MatrixXd>& dataA,
const QVector<MatrixXd>& dataB)
138 const int nSubjects = dataA.size();
139 const int nChannels =
static_cast<int>(dataA[0].rows());
140 const int nTimes =
static_cast<int>(dataA[0].cols());
144 MatrixXd sumDiff = MatrixXd::Zero(nChannels, nTimes);
145 MatrixXd sumDiffSq = MatrixXd::Zero(nChannels, nTimes);
147 for (
int s = 0; s < nSubjects; ++s) {
148 MatrixXd diff = dataA[s] - dataB[s];
150 sumDiffSq += diff.cwiseProduct(diff);
153 double n =
static_cast<double>(nSubjects);
154 MatrixXd mean = sumDiff / n;
155 MatrixXd variance = (sumDiffSq - sumDiff.cwiseProduct(sumDiff) / n) / (n - 1.0);
158 variance = variance.cwiseMax(1.0e-30);
160 MatrixXd tMap = mean.array() / (variance.array().sqrt() / std::sqrt(n));
166QPair<MatrixXi, QVector<double>> StatsCluster::findClusters(
167 const MatrixXd& tMap,
169 const SparseMatrix<int>& adjacency,
172 const int nChannels =
static_cast<int>(tMap.rows());
173 const int nTimes =
static_cast<int>(tMap.cols());
175 MatrixXi clusterIds = MatrixXi::Zero(nChannels, nTimes);
176 QVector<double> clusterStats;
177 int currentClusterId = 0;
180 auto bfsClusters = [&](
bool positive) {
182 MatrixXi visited = MatrixXi::Zero(nChannels, nTimes);
184 for (
int ch = 0; ch < nChannels; ++ch) {
185 for (
int t = 0; t < nTimes; ++t) {
186 double val = tMap(ch, t);
187 bool suprathreshold = positive ? (val > threshold) : (val < -threshold);
189 if (!suprathreshold || visited(ch, t))
continue;
193 double clusterSum = 0.0;
194 std::queue<std::pair<int, int>> queue;
198 while (!queue.empty()) {
199 auto [curCh, curT] = queue.front();
202 clusterIds(curCh, curT) = positive ? currentClusterId : -currentClusterId;
203 clusterSum += tMap(curCh, curT);
206 for (
int dt = -1; dt <= 1; dt += 2) {
208 if (nt < 0 || nt >= nTimes)
continue;
209 double nval = tMap(curCh, nt);
210 bool nSupra = positive ? (nval > threshold) : (nval < -threshold);
211 if (nSupra && !visited(curCh, nt)) {
212 visited(curCh, nt) = 1;
213 queue.push({curCh, nt});
218 for (SparseMatrix<int>::InnerIterator it(adjacency, curCh); it; ++it) {
219 int nCh =
static_cast<int>(it.row());
220 double nval = tMap(nCh, curT);
221 bool nSupra = positive ? (nval > threshold) : (nval < -threshold);
222 if (nSupra && !visited(nCh, curT)) {
223 visited(nCh, curT) = 1;
224 queue.push({nCh, curT});
229 clusterStats.append(clusterSum);
241 return {clusterIds, clusterStats};
246double StatsCluster::permuteOnce(
247 const QVector<MatrixXd>& allData,
249 const SparseMatrix<int>& adjacency,
253 const int nTotal = allData.size();
256 std::vector<int> indices(nTotal);
257 std::iota(indices.begin(), indices.end(), 0);
260 QRandomGenerator rng(QRandomGenerator::global()->generate());
261 for (
int i = nTotal - 1; i > 0; --i) {
262 int j =
static_cast<int>(rng.bounded(i + 1));
263 std::swap(indices[i], indices[j]);
267 QVector<MatrixXd> permA, permB;
269 permB.reserve(nTotal - nA);
270 for (
int i = 0; i < nTotal; ++i) {
272 permA.append(allData[indices[i]]);
274 permB.append(allData[indices[i]]);
279 MatrixXd tMap = computeTMap(permA, permB);
280 auto [clusterIds, clusterStats] = findClusters(tMap, threshold, adjacency, tail);
283 double maxStat = 0.0;
284 for (
double s : clusterStats) {
285 double absS = std::fabs(s);
286 if (absS > maxStat) {
295double StatsCluster::inverseTCdf(
double p,
int df)
302 double targetCdf = 1.0 - p / 2.0;
307 for (
int iter = 0; iter < 100; ++iter) {
308 double mid = (lo + hi) / 2.0;
310 if (cdf < targetCdf) {
317 return (lo + hi) / 2.0;
322MatrixXd StatsCluster::computeOneSampleTMap(
const QVector<MatrixXd>& data)
324 const int nSubjects = data.size();
325 const int nVertices =
static_cast<int>(data[0].rows());
326 const int nTimes =
static_cast<int>(data[0].cols());
328 MatrixXd sum = MatrixXd::Zero(nVertices, nTimes);
329 MatrixXd sumSq = MatrixXd::Zero(nVertices, nTimes);
331 for (
int s = 0; s < nSubjects; ++s) {
333 sumSq += data[s].cwiseProduct(data[s]);
336 double n =
static_cast<double>(nSubjects);
337 MatrixXd mean = sum / n;
338 MatrixXd variance = (sumSq - sum.cwiseProduct(sum) / n) / (n - 1.0);
339 variance = variance.cwiseMax(1.0e-30);
341 MatrixXd tMap = mean.array() / (variance.array().sqrt() / std::sqrt(n));
347MatrixXd StatsCluster::computeFMap(
const QVector<QVector<MatrixXd>>& conditions)
349 const int nConditions = conditions.size();
350 const int nVertices =
static_cast<int>(conditions[0][0].rows());
351 const int nTimes =
static_cast<int>(conditions[0][0].cols());
355 for (
int c = 0; c < nConditions; ++c) {
356 nTotal += conditions[c].size();
360 MatrixXd grandSum = MatrixXd::Zero(nVertices, nTimes);
361 for (
int c = 0; c < nConditions; ++c) {
362 for (
int s = 0; s < conditions[c].size(); ++s) {
363 grandSum += conditions[c][s];
366 MatrixXd grandMean = grandSum /
static_cast<double>(nTotal);
369 MatrixXd ssBetween = MatrixXd::Zero(nVertices, nTimes);
370 MatrixXd ssWithin = MatrixXd::Zero(nVertices, nTimes);
372 for (
int c = 0; c < nConditions; ++c) {
373 int nc = conditions[c].size();
374 MatrixXd condSum = MatrixXd::Zero(nVertices, nTimes);
375 for (
int s = 0; s < nc; ++s) {
376 condSum += conditions[c][s];
378 MatrixXd condMean = condSum /
static_cast<double>(nc);
380 MatrixXd diff = condMean - grandMean;
381 ssBetween +=
static_cast<double>(nc) * diff.cwiseProduct(diff);
383 for (
int s = 0; s < nc; ++s) {
384 MatrixXd residual = conditions[c][s] - condMean;
385 ssWithin += residual.cwiseProduct(residual);
389 int dfBetween = nConditions - 1;
390 int dfWithin = nTotal - nConditions;
392 MatrixXd msBetween = ssBetween /
static_cast<double>(dfBetween);
393 MatrixXd msWithin = ssWithin /
static_cast<double>(dfWithin);
394 msWithin = msWithin.cwiseMax(1.0e-30);
396 MatrixXd fMap = msBetween.array() / msWithin.array();
402QPair<MatrixXi, QVector<double>> StatsCluster::findClustersFlat(
403 const MatrixXd& statMap,
405 const SparseMatrix<int>& adjacency,
408 const int nVertices =
static_cast<int>(statMap.rows());
409 const int nTimes =
static_cast<int>(statMap.cols());
410 const int nTotal = nVertices * nTimes;
412 MatrixXi clusterIds = MatrixXi::Zero(nVertices, nTimes);
413 QVector<double> clusterStats;
414 int currentClusterId = 0;
417 std::vector<bool> visited(nTotal,
false);
419 for (
int idx = 0; idx < nTotal; ++idx) {
420 int v = idx / nTimes;
421 int t = idx % nTimes;
422 double val = statMap(v, t);
424 bool suprathreshold = positiveOnly ? (val > threshold) : (val < -threshold);
425 if (!suprathreshold || visited[idx])
continue;
428 double clusterSum = 0.0;
429 std::queue<int> queue;
433 while (!queue.empty()) {
434 int curIdx = queue.front();
437 int curV = curIdx / nTimes;
438 int curT = curIdx % nTimes;
440 clusterIds(curV, curT) = positiveOnly ? currentClusterId : -currentClusterId;
441 clusterSum += statMap(curV, curT);
444 for (SparseMatrix<int>::InnerIterator it(adjacency, curIdx); it; ++it) {
445 int nIdx = static_cast<int>(it.row());
446 if (visited[nIdx]) continue;
447 int nV = nIdx / nTimes;
448 int nT = nIdx % nTimes;
449 double nval = statMap(nV, nT);
450 bool nSupra = positiveOnly ? (nval > threshold) : (nval < -threshold);
452 visited[nIdx] = true;
458 clusterStats.append(clusterSum);
461 return {clusterIds, clusterStats};
466double StatsCluster::permuteOnceOneSample(
467 const QVector<MatrixXd>& data,
468 const SparseMatrix<int>& adjacency,
472 const int nSubjects = data.size();
475 QRandomGenerator rng(QRandomGenerator::global()->generate());
476 QVector<int> signs(nSubjects);
477 for (
int s = 0; s < nSubjects; ++s) {
478 signs[s] = rng.bounded(2) == 0 ? 1 : -1;
482 QVector<MatrixXd> flipped(nSubjects);
483 for (
int s = 0; s < nSubjects; ++s) {
484 flipped[s] = data[s] *
static_cast<double>(signs[s]);
488 MatrixXd tMap = computeOneSampleTMap(flipped);
490 double maxStat = 0.0;
493 auto [ids, stats] = findClustersFlat(tMap, threshold, adjacency,
true);
494 for (
double s : stats) {
495 if (std::fabs(s) > maxStat) maxStat = std::fabs(s);
499 auto [ids, stats] = findClustersFlat(tMap, threshold, adjacency,
false);
500 for (
double s : stats) {
501 if (std::fabs(s) > maxStat) maxStat = std::fabs(s);
510double StatsCluster::permuteOnceFTest(
511 const QVector<MatrixXd>& allData,
512 const QVector<int>& groupSizes,
513 const SparseMatrix<int>& adjacency,
516 const int nTotal = allData.size();
519 std::vector<int> indices(nTotal);
520 std::iota(indices.begin(), indices.end(), 0);
522 QRandomGenerator rng(QRandomGenerator::global()->generate());
523 for (
int i = nTotal - 1; i > 0; --i) {
524 int j =
static_cast<int>(rng.bounded(i + 1));
525 std::swap(indices[i], indices[j]);
529 QVector<QVector<MatrixXd>> permConditions;
531 for (
int c = 0; c < groupSizes.size(); ++c) {
532 QVector<MatrixXd> group;
533 group.reserve(groupSizes[c]);
534 for (
int s = 0; s < groupSizes[c]; ++s) {
535 group.append(allData[indices[offset + s]]);
537 permConditions.append(group);
538 offset += groupSizes[c];
542 MatrixXd fMap = computeFMap(permConditions);
543 auto [ids, stats] = findClustersFlat(fMap, threshold, adjacency,
true);
545 double maxStat = 0.0;
546 for (
double s : stats) {
547 if (s > maxStat) maxStat = s;
555 const QVector<MatrixXd>& data,
556 const SparseMatrix<int>& adjacency,
562 MatrixXd tObs = computeOneSampleTMap(data);
565 MatrixXi clusterIdsPos, clusterIdsNeg;
566 QVector<double> clusterStats;
569 auto [ids, stats] = findClustersFlat(tObs, threshold, adjacency,
true);
571 clusterStats.append(stats);
574 auto [ids, stats] = findClustersFlat(tObs, threshold, adjacency,
false);
576 clusterStats.append(stats);
580 MatrixXi clusterIds = MatrixXi::Zero(tObs.rows(), tObs.cols());
581 if (clusterIdsPos.size() > 0) clusterIds += clusterIdsPos;
582 if (clusterIdsNeg.size() > 0) clusterIds += clusterIdsNeg;
585 QVector<int> permIndices(nPermutations);
586 std::iota(permIndices.begin(), permIndices.end(), 0);
588 std::function<double(
int)> permuteFunc = [&](int) ->
double {
589 return permuteOnceOneSample(data, adjacency, threshold, tail);
592 QFuture<double> future = QtConcurrent::mapped(permIndices, permuteFunc);
593 future.waitForFinished();
595 QVector<double> nullDist(nPermutations);
596 for (
int i = 0; i < nPermutations; ++i) {
597 nullDist[i] = future.resultAt(i);
599 std::sort(nullDist.begin(), nullDist.end());
602 QVector<double> clusterPvals(clusterStats.size());
603 for (
int c = 0; c < clusterStats.size(); ++c) {
604 double obsStat = std::fabs(clusterStats[c]);
606 for (
int p = 0; p < nPermutations; ++p) {
607 if (nullDist[p] >= obsStat) {
611 clusterPvals[c] =
static_cast<double>(count) /
static_cast<double>(nPermutations);
626 const QVector<QVector<MatrixXd>>& conditions,
627 const SparseMatrix<int>& adjacency,
632 MatrixXd fObs = computeFMap(conditions);
635 auto [clusterIds, clusterStats] = findClustersFlat(fObs, threshold, adjacency,
true);
638 QVector<MatrixXd> allData;
639 QVector<int> groupSizes;
640 for (
int c = 0; c < conditions.size(); ++c) {
641 groupSizes.append(conditions[c].size());
642 for (
int s = 0; s < conditions[c].size(); ++s) {
643 allData.append(conditions[c][s]);
648 QVector<int> permIndices(nPermutations);
649 std::iota(permIndices.begin(), permIndices.end(), 0);
651 std::function<double(
int)> permuteFunc = [&](int) ->
double {
652 return permuteOnceFTest(allData, groupSizes, adjacency, threshold);
655 QFuture<double> future = QtConcurrent::mapped(permIndices, permuteFunc);
656 future.waitForFinished();
658 QVector<double> nullDist(nPermutations);
659 for (
int i = 0; i < nPermutations; ++i) {
660 nullDist[i] = future.resultAt(i);
662 std::sort(nullDist.begin(), nullDist.end());
665 QVector<double> clusterPvals(clusterStats.size());
666 for (
int c = 0; c < clusterStats.size(); ++c) {
667 double obsStat = clusterStats[c];
669 for (
int p = 0; p < nPermutations; ++p) {
670 if (nullDist[p] >= obsStat) {
674 clusterPvals[c] =
static_cast<double>(count) /
static_cast<double>(nPermutations);
689 const MatrixXd& statMap,
690 const SparseMatrix<int>& adjacency,
695 const int nVertices =
static_cast<int>(statMap.rows());
696 const int nTimes =
static_cast<int>(statMap.cols());
697 const int nTotal = nVertices * nTimes;
699 MatrixXd tfceMap = MatrixXd::Zero(nVertices, nTimes);
702 auto runTfce = [&](
const MatrixXd& absMap,
bool isPositive) {
703 double maxVal = absMap.maxCoeff();
704 if (maxVal <= 0.0)
return;
706 double dh = maxVal /
static_cast<double>(nSteps);
708 for (
int step = 1; step <= nSteps; ++step) {
709 double h = dh *
static_cast<double>(step);
712 std::vector<bool> visited(nTotal,
false);
714 for (
int idx = 0; idx < nTotal; ++idx) {
715 int v = idx / nTimes;
716 int t = idx % nTimes;
718 if (absMap(v, t) < h || visited[idx])
continue;
721 std::vector<int> cluster;
722 std::queue<int> queue;
726 while (!queue.empty()) {
727 int curIdx = queue.front();
729 cluster.push_back(curIdx);
731 for (SparseMatrix<int>::InnerIterator it(adjacency, curIdx); it; ++it) {
732 int nIdx =
static_cast<int>(it.row());
733 if (visited[nIdx])
continue;
734 int nV = nIdx / nTimes;
735 int nT = nIdx % nTimes;
736 if (absMap(nV, nT) >= h) {
737 visited[nIdx] =
true;
744 double extent =
static_cast<double>(cluster.size());
745 double contribution = std::pow(extent, E) * std::pow(h, H) * dh;
747 for (
int cIdx : cluster) {
748 int cv = cIdx / nTimes;
749 int ct = cIdx % nTimes;
751 tfceMap(cv, ct) += contribution;
753 tfceMap(cv, ct) -= contribution;
761 MatrixXd posMap = statMap.cwiseMax(0.0);
762 runTfce(posMap,
true);
765 MatrixXd negMap = (-statMap).cwiseMax(0.0);
766 runTfce(negMap,
false);
StatsCluster class declaration.
StatsTtest class declaration.
StatsFtest class declaration.
Statistical testing (t-tests, F-tests, cluster permutation, multiple comparison correction).
Eigen::MatrixXi matClusterIds
QVector< double > vecClusterStats
QVector< double > vecClusterPvals
static StatsClusterResult permutationTest(const QVector< Eigen::MatrixXd > &dataA, const QVector< Eigen::MatrixXd > &dataB, const Eigen::SparseMatrix< int > &adjacency, int nPermutations=1024, double clusterAlpha=0.05, double pThreshold=0.05, StatsTailType tail=StatsTailType::Both)
static Eigen::MatrixXd tfce(const Eigen::MatrixXd &statMap, const Eigen::SparseMatrix< int > &adjacency, double E=0.5, double H=2.0, int nSteps=100)
static StatsClusterResult oneSamplePermutationTest(const QVector< Eigen::MatrixXd > &data, const Eigen::SparseMatrix< int > &adjacency, double threshold, int nPermutations, StatsTailType tail)
static StatsClusterResult fTestPermutationTest(const QVector< QVector< Eigen::MatrixXd > > &conditions, const Eigen::SparseMatrix< int > &adjacency, double threshold, int nPermutations)
static double tCdf(double t, int df)