v2.0.0
Loading...
Searching...
No Matches
sts_cluster.cpp
Go to the documentation of this file.
1//=============================================================================================================
34
35//=============================================================================================================
36// INCLUDES
37//=============================================================================================================
38
39#include "sts_cluster.h"
40#include "sts_ttest.h"
41#include "sts_ftest.h"
42
43#include <QtConcurrent>
44#include <QRandomGenerator>
45
46#include <cmath>
47#include <algorithm>
48#include <queue>
49#include <vector>
50
51//=============================================================================================================
52// USED NAMESPACES
53//=============================================================================================================
54
55using namespace STSLIB;
56using namespace Eigen;
57
58//=============================================================================================================
59// DEFINE METHODS
60//=============================================================================================================
61
63 const QVector<MatrixXd>& dataA,
64 const QVector<MatrixXd>& dataB,
65 const SparseMatrix<int>& adjacency,
66 int nPermutations,
67 double clusterAlpha,
68 double pThreshold,
69 StatsTailType tail)
70{
71 Q_UNUSED(pThreshold);
72
73 const int nA = dataA.size();
74
75 // Step 1: Compute the observed t-map
76 MatrixXd tObs = computeTMap(dataA, dataB);
77
78 // Step 2: Determine cluster-forming threshold from clusterAlpha and df
79 int df = nA - 1;
80 double threshold = inverseTCdf(clusterAlpha, df);
81
82 // Step 3: Find observed clusters
83 auto [clusterIds, clusterStats] = findClusters(tObs, threshold, adjacency, tail);
84
85 // Step 4: Combine all data for permutation
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);
90
91 // Step 5: Build null distribution via permutations (using QtConcurrent)
92 QVector<double> nullDist(nPermutations);
93
94 // Perform permutations in parallel
95 QVector<int> permIndices(nPermutations);
96 std::iota(permIndices.begin(), permIndices.end(), 0);
97
98 std::function<double(int)> permuteFunc = [&](int) -> double {
99 return permuteOnce(allData, nA, adjacency, threshold, tail);
100 };
101
102 QFuture<double> future = QtConcurrent::mapped(permIndices, permuteFunc);
103 future.waitForFinished();
104
105 for (int i = 0; i < nPermutations; ++i) {
106 nullDist[i] = future.resultAt(i);
107 }
108
109 // Sort null distribution for percentile computation
110 std::sort(nullDist.begin(), nullDist.end());
111
112 // Step 6: Compute cluster p-values
113 QVector<double> clusterPvals(clusterStats.size());
114 for (int c = 0; c < clusterStats.size(); ++c) {
115 double obsStat = std::fabs(clusterStats[c]);
116 int count = 0;
117 for (int p = 0; p < nPermutations; ++p) {
118 if (nullDist[p] >= obsStat) {
119 count++;
120 }
121 }
122 clusterPvals[c] = static_cast<double>(count) / static_cast<double>(nPermutations);
123 }
124
125 StatsClusterResult result;
126 result.matTObs = tObs;
127 result.vecClusterStats = clusterStats;
128 result.vecClusterPvals = clusterPvals;
129 result.matClusterIds = clusterIds;
130 result.clusterThreshold = threshold;
131 return result;
132}
133
134//=============================================================================================================
135
136MatrixXd StatsCluster::computeTMap(const QVector<MatrixXd>& dataA, const QVector<MatrixXd>& dataB)
137{
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());
141
142 // Compute difference for each subject
143 // Then compute paired t-statistic at each (channel, time)
144 MatrixXd sumDiff = MatrixXd::Zero(nChannels, nTimes);
145 MatrixXd sumDiffSq = MatrixXd::Zero(nChannels, nTimes);
146
147 for (int s = 0; s < nSubjects; ++s) {
148 MatrixXd diff = dataA[s] - dataB[s];
149 sumDiff += diff;
150 sumDiffSq += diff.cwiseProduct(diff);
151 }
152
153 double n = static_cast<double>(nSubjects);
154 MatrixXd mean = sumDiff / n;
155 MatrixXd variance = (sumDiffSq - sumDiff.cwiseProduct(sumDiff) / n) / (n - 1.0);
156
157 // Clamp variance to avoid division by zero
158 variance = variance.cwiseMax(1.0e-30);
159
160 MatrixXd tMap = mean.array() / (variance.array().sqrt() / std::sqrt(n));
161 return tMap;
162}
163
164//=============================================================================================================
165
166QPair<MatrixXi, QVector<double>> StatsCluster::findClusters(
167 const MatrixXd& tMap,
168 double threshold,
169 const SparseMatrix<int>& adjacency,
170 StatsTailType tail)
171{
172 const int nChannels = static_cast<int>(tMap.rows());
173 const int nTimes = static_cast<int>(tMap.cols());
174
175 MatrixXi clusterIds = MatrixXi::Zero(nChannels, nTimes);
176 QVector<double> clusterStats;
177 int currentClusterId = 0;
178
179 // Helper lambda: BFS on combined spatial+temporal adjacency for one polarity
180 auto bfsClusters = [&](bool positive) {
181 // Build a visited mask
182 MatrixXi visited = MatrixXi::Zero(nChannels, nTimes);
183
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);
188
189 if (!suprathreshold || visited(ch, t)) continue;
190
191 // Start BFS for a new cluster
192 currentClusterId++;
193 double clusterSum = 0.0;
194 std::queue<std::pair<int, int>> queue;
195 queue.push({ch, t});
196 visited(ch, t) = 1;
197
198 while (!queue.empty()) {
199 auto [curCh, curT] = queue.front();
200 queue.pop();
201
202 clusterIds(curCh, curT) = positive ? currentClusterId : -currentClusterId;
203 clusterSum += tMap(curCh, curT);
204
205 // Temporal neighbors (consecutive time points)
206 for (int dt = -1; dt <= 1; dt += 2) {
207 int nt = curT + dt;
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});
214 }
215 }
216
217 // Spatial neighbors at the same time point
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});
225 }
226 }
227 }
228
229 clusterStats.append(clusterSum);
230 }
231 }
232 };
233
234 if (tail == StatsTailType::Both || tail == StatsTailType::Right) {
235 bfsClusters(true);
236 }
237 if (tail == StatsTailType::Both || tail == StatsTailType::Left) {
238 bfsClusters(false);
239 }
240
241 return {clusterIds, clusterStats};
242}
243
244//=============================================================================================================
245
246double StatsCluster::permuteOnce(
247 const QVector<MatrixXd>& allData,
248 int nA,
249 const SparseMatrix<int>& adjacency,
250 double threshold,
251 StatsTailType tail)
252{
253 const int nTotal = allData.size();
254
255 // Generate a random permutation of indices
256 std::vector<int> indices(nTotal);
257 std::iota(indices.begin(), indices.end(), 0);
258
259 // Use thread-safe random generator
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]);
264 }
265
266 // Split into permuted groups
267 QVector<MatrixXd> permA, permB;
268 permA.reserve(nA);
269 permB.reserve(nTotal - nA);
270 for (int i = 0; i < nTotal; ++i) {
271 if (i < nA) {
272 permA.append(allData[indices[i]]);
273 } else {
274 permB.append(allData[indices[i]]);
275 }
276 }
277
278 // Compute t-map and find clusters
279 MatrixXd tMap = computeTMap(permA, permB);
280 auto [clusterIds, clusterStats] = findClusters(tMap, threshold, adjacency, tail);
281
282 // Return the maximum absolute cluster statistic
283 double maxStat = 0.0;
284 for (double s : clusterStats) {
285 double absS = std::fabs(s);
286 if (absS > maxStat) {
287 maxStat = absS;
288 }
289 }
290 return maxStat;
291}
292
293//=============================================================================================================
294
295double StatsCluster::inverseTCdf(double p, int df)
296{
297 // Approximate inverse t CDF for two-tailed threshold.
298 // We want the t-value such that P(|T| > t) = p, i.e., P(T > t) = p/2.
299 // This means we want the (1 - p/2) quantile.
300
301 // Use a simple bisection on the tCdf function.
302 double targetCdf = 1.0 - p / 2.0;
303
304 double lo = 0.0;
305 double hi = 100.0;
306
307 for (int iter = 0; iter < 100; ++iter) {
308 double mid = (lo + hi) / 2.0;
309 double cdf = StatsTtest::tCdf(mid, df);
310 if (cdf < targetCdf) {
311 lo = mid;
312 } else {
313 hi = mid;
314 }
315 }
316
317 return (lo + hi) / 2.0;
318}
319
320//=============================================================================================================
321
322MatrixXd StatsCluster::computeOneSampleTMap(const QVector<MatrixXd>& data)
323{
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());
327
328 MatrixXd sum = MatrixXd::Zero(nVertices, nTimes);
329 MatrixXd sumSq = MatrixXd::Zero(nVertices, nTimes);
330
331 for (int s = 0; s < nSubjects; ++s) {
332 sum += data[s];
333 sumSq += data[s].cwiseProduct(data[s]);
334 }
335
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);
340
341 MatrixXd tMap = mean.array() / (variance.array().sqrt() / std::sqrt(n));
342 return tMap;
343}
344
345//=============================================================================================================
346
347MatrixXd StatsCluster::computeFMap(const QVector<QVector<MatrixXd>>& conditions)
348{
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());
352
353 // Count total subjects and build per-condition means
354 int nTotal = 0;
355 for (int c = 0; c < nConditions; ++c) {
356 nTotal += conditions[c].size();
357 }
358
359 // Grand mean
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];
364 }
365 }
366 MatrixXd grandMean = grandSum / static_cast<double>(nTotal);
367
368 // SS between and SS within
369 MatrixXd ssBetween = MatrixXd::Zero(nVertices, nTimes);
370 MatrixXd ssWithin = MatrixXd::Zero(nVertices, nTimes);
371
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];
377 }
378 MatrixXd condMean = condSum / static_cast<double>(nc);
379
380 MatrixXd diff = condMean - grandMean;
381 ssBetween += static_cast<double>(nc) * diff.cwiseProduct(diff);
382
383 for (int s = 0; s < nc; ++s) {
384 MatrixXd residual = conditions[c][s] - condMean;
385 ssWithin += residual.cwiseProduct(residual);
386 }
387 }
388
389 int dfBetween = nConditions - 1;
390 int dfWithin = nTotal - nConditions;
391
392 MatrixXd msBetween = ssBetween / static_cast<double>(dfBetween);
393 MatrixXd msWithin = ssWithin / static_cast<double>(dfWithin);
394 msWithin = msWithin.cwiseMax(1.0e-30);
395
396 MatrixXd fMap = msBetween.array() / msWithin.array();
397 return fMap;
398}
399
400//=============================================================================================================
401
402QPair<MatrixXi, QVector<double>> StatsCluster::findClustersFlat(
403 const MatrixXd& statMap,
404 double threshold,
405 const SparseMatrix<int>& adjacency,
406 bool positiveOnly)
407{
408 const int nVertices = static_cast<int>(statMap.rows());
409 const int nTimes = static_cast<int>(statMap.cols());
410 const int nTotal = nVertices * nTimes;
411
412 MatrixXi clusterIds = MatrixXi::Zero(nVertices, nTimes);
413 QVector<double> clusterStats;
414 int currentClusterId = 0;
415
416 // BFS on the flat spatio-temporal adjacency
417 std::vector<bool> visited(nTotal, false);
418
419 for (int idx = 0; idx < nTotal; ++idx) {
420 int v = idx / nTimes;
421 int t = idx % nTimes;
422 double val = statMap(v, t);
423
424 bool suprathreshold = positiveOnly ? (val > threshold) : (val < -threshold);
425 if (!suprathreshold || visited[idx]) continue;
426
427 currentClusterId++;
428 double clusterSum = 0.0;
429 std::queue<int> queue;
430 queue.push(idx);
431 visited[idx] = true;
432
433 while (!queue.empty()) {
434 int curIdx = queue.front();
435 queue.pop();
436
437 int curV = curIdx / nTimes;
438 int curT = curIdx % nTimes;
439
440 clusterIds(curV, curT) = positiveOnly ? currentClusterId : -currentClusterId;
441 clusterSum += statMap(curV, curT);
442
443 // Iterate over adjacency neighbors
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);
451 if (nSupra) {
452 visited[nIdx] = true;
453 queue.push(nIdx);
454 }
455 }
456 }
457
458 clusterStats.append(clusterSum);
459 }
460
461 return {clusterIds, clusterStats};
462}
463
464//=============================================================================================================
465
466double StatsCluster::permuteOnceOneSample(
467 const QVector<MatrixXd>& data,
468 const SparseMatrix<int>& adjacency,
469 double threshold,
470 StatsTailType tail)
471{
472 const int nSubjects = data.size();
473
474 // Generate random sign flips
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;
479 }
480
481 // Apply sign flips
482 QVector<MatrixXd> flipped(nSubjects);
483 for (int s = 0; s < nSubjects; ++s) {
484 flipped[s] = data[s] * static_cast<double>(signs[s]);
485 }
486
487 // Compute t-map and find clusters
488 MatrixXd tMap = computeOneSampleTMap(flipped);
489
490 double maxStat = 0.0;
491
492 if (tail == StatsTailType::Both || tail == StatsTailType::Right) {
493 auto [ids, stats] = findClustersFlat(tMap, threshold, adjacency, true);
494 for (double s : stats) {
495 if (std::fabs(s) > maxStat) maxStat = std::fabs(s);
496 }
497 }
498 if (tail == StatsTailType::Both || tail == StatsTailType::Left) {
499 auto [ids, stats] = findClustersFlat(tMap, threshold, adjacency, false);
500 for (double s : stats) {
501 if (std::fabs(s) > maxStat) maxStat = std::fabs(s);
502 }
503 }
504
505 return maxStat;
506}
507
508//=============================================================================================================
509
510double StatsCluster::permuteOnceFTest(
511 const QVector<MatrixXd>& allData,
512 const QVector<int>& groupSizes,
513 const SparseMatrix<int>& adjacency,
514 double threshold)
515{
516 const int nTotal = allData.size();
517
518 // Random permutation of indices
519 std::vector<int> indices(nTotal);
520 std::iota(indices.begin(), indices.end(), 0);
521
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]);
526 }
527
528 // Rebuild condition groups
529 QVector<QVector<MatrixXd>> permConditions;
530 int offset = 0;
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]]);
536 }
537 permConditions.append(group);
538 offset += groupSizes[c];
539 }
540
541 // Compute F-map and find clusters (F is always positive)
542 MatrixXd fMap = computeFMap(permConditions);
543 auto [ids, stats] = findClustersFlat(fMap, threshold, adjacency, true);
544
545 double maxStat = 0.0;
546 for (double s : stats) {
547 if (s > maxStat) maxStat = s;
548 }
549 return maxStat;
550}
551
552//=============================================================================================================
553
555 const QVector<MatrixXd>& data,
556 const SparseMatrix<int>& adjacency,
557 double threshold,
558 int nPermutations,
559 StatsTailType tail)
560{
561 // Step 1: Compute observed one-sample t-map
562 MatrixXd tObs = computeOneSampleTMap(data);
563
564 // Step 2: Find observed clusters
565 MatrixXi clusterIdsPos, clusterIdsNeg;
566 QVector<double> clusterStats;
567
568 if (tail == StatsTailType::Both || tail == StatsTailType::Right) {
569 auto [ids, stats] = findClustersFlat(tObs, threshold, adjacency, true);
570 clusterIdsPos = ids;
571 clusterStats.append(stats);
572 }
573 if (tail == StatsTailType::Both || tail == StatsTailType::Left) {
574 auto [ids, stats] = findClustersFlat(tObs, threshold, adjacency, false);
575 clusterIdsNeg = ids;
576 clusterStats.append(stats);
577 }
578
579 // Merge cluster ID maps
580 MatrixXi clusterIds = MatrixXi::Zero(tObs.rows(), tObs.cols());
581 if (clusterIdsPos.size() > 0) clusterIds += clusterIdsPos;
582 if (clusterIdsNeg.size() > 0) clusterIds += clusterIdsNeg;
583
584 // Step 3: Build null distribution via sign-flip permutations
585 QVector<int> permIndices(nPermutations);
586 std::iota(permIndices.begin(), permIndices.end(), 0);
587
588 std::function<double(int)> permuteFunc = [&](int) -> double {
589 return permuteOnceOneSample(data, adjacency, threshold, tail);
590 };
591
592 QFuture<double> future = QtConcurrent::mapped(permIndices, permuteFunc);
593 future.waitForFinished();
594
595 QVector<double> nullDist(nPermutations);
596 for (int i = 0; i < nPermutations; ++i) {
597 nullDist[i] = future.resultAt(i);
598 }
599 std::sort(nullDist.begin(), nullDist.end());
600
601 // Step 4: Compute cluster p-values
602 QVector<double> clusterPvals(clusterStats.size());
603 for (int c = 0; c < clusterStats.size(); ++c) {
604 double obsStat = std::fabs(clusterStats[c]);
605 int count = 0;
606 for (int p = 0; p < nPermutations; ++p) {
607 if (nullDist[p] >= obsStat) {
608 count++;
609 }
610 }
611 clusterPvals[c] = static_cast<double>(count) / static_cast<double>(nPermutations);
612 }
613
614 StatsClusterResult result;
615 result.matTObs = tObs;
616 result.vecClusterStats = clusterStats;
617 result.vecClusterPvals = clusterPvals;
618 result.matClusterIds = clusterIds;
619 result.clusterThreshold = threshold;
620 return result;
621}
622
623//=============================================================================================================
624
626 const QVector<QVector<MatrixXd>>& conditions,
627 const SparseMatrix<int>& adjacency,
628 double threshold,
629 int nPermutations)
630{
631 // Step 1: Compute observed F-map
632 MatrixXd fObs = computeFMap(conditions);
633
634 // Step 2: Find observed clusters (F is always positive, one-tailed)
635 auto [clusterIds, clusterStats] = findClustersFlat(fObs, threshold, adjacency, true);
636
637 // Step 3: Flatten all data and record group sizes
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]);
644 }
645 }
646
647 // Step 4: Build null distribution via label-shuffle permutations
648 QVector<int> permIndices(nPermutations);
649 std::iota(permIndices.begin(), permIndices.end(), 0);
650
651 std::function<double(int)> permuteFunc = [&](int) -> double {
652 return permuteOnceFTest(allData, groupSizes, adjacency, threshold);
653 };
654
655 QFuture<double> future = QtConcurrent::mapped(permIndices, permuteFunc);
656 future.waitForFinished();
657
658 QVector<double> nullDist(nPermutations);
659 for (int i = 0; i < nPermutations; ++i) {
660 nullDist[i] = future.resultAt(i);
661 }
662 std::sort(nullDist.begin(), nullDist.end());
663
664 // Step 5: Compute cluster p-values
665 QVector<double> clusterPvals(clusterStats.size());
666 for (int c = 0; c < clusterStats.size(); ++c) {
667 double obsStat = clusterStats[c];
668 int count = 0;
669 for (int p = 0; p < nPermutations; ++p) {
670 if (nullDist[p] >= obsStat) {
671 count++;
672 }
673 }
674 clusterPvals[c] = static_cast<double>(count) / static_cast<double>(nPermutations);
675 }
676
677 StatsClusterResult result;
678 result.matTObs = fObs;
679 result.vecClusterStats = clusterStats;
680 result.vecClusterPvals = clusterPvals;
681 result.matClusterIds = clusterIds;
682 result.clusterThreshold = threshold;
683 return result;
684}
685
686//=============================================================================================================
687
689 const MatrixXd& statMap,
690 const SparseMatrix<int>& adjacency,
691 double E,
692 double H,
693 int nSteps)
694{
695 const int nVertices = static_cast<int>(statMap.rows());
696 const int nTimes = static_cast<int>(statMap.cols());
697 const int nTotal = nVertices * nTimes;
698
699 MatrixXd tfceMap = MatrixXd::Zero(nVertices, nTimes);
700
701 // Helper: run TFCE on one polarity
702 auto runTfce = [&](const MatrixXd& absMap, bool isPositive) {
703 double maxVal = absMap.maxCoeff();
704 if (maxVal <= 0.0) return;
705
706 double dh = maxVal / static_cast<double>(nSteps);
707
708 for (int step = 1; step <= nSteps; ++step) {
709 double h = dh * static_cast<double>(step);
710
711 // Find connected components above threshold h
712 std::vector<bool> visited(nTotal, false);
713
714 for (int idx = 0; idx < nTotal; ++idx) {
715 int v = idx / nTimes;
716 int t = idx % nTimes;
717
718 if (absMap(v, t) < h || visited[idx]) continue;
719
720 // BFS to find cluster
721 std::vector<int> cluster;
722 std::queue<int> queue;
723 queue.push(idx);
724 visited[idx] = true;
725
726 while (!queue.empty()) {
727 int curIdx = queue.front();
728 queue.pop();
729 cluster.push_back(curIdx);
730
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;
738 queue.push(nIdx);
739 }
740 }
741 }
742
743 // Compute contribution: e^E * h^H * dh
744 double extent = static_cast<double>(cluster.size());
745 double contribution = std::pow(extent, E) * std::pow(h, H) * dh;
746
747 for (int cIdx : cluster) {
748 int cv = cIdx / nTimes;
749 int ct = cIdx % nTimes;
750 if (isPositive) {
751 tfceMap(cv, ct) += contribution;
752 } else {
753 tfceMap(cv, ct) -= contribution;
754 }
755 }
756 }
757 }
758 };
759
760 // Positive values
761 MatrixXd posMap = statMap.cwiseMax(0.0);
762 runTfce(posMap, true);
763
764 // Negative values (use absolute values, then negate contributions)
765 MatrixXd negMap = (-statMap).cwiseMax(0.0);
766 runTfce(negMap, false);
767
768 return tfceMap;
769}
StatsCluster class declaration.
StatsTtest class declaration.
StatsFtest class declaration.
Statistical testing (t-tests, F-tests, cluster permutation, multiple comparison correction).
StatsTailType
Definition sts_types.h:49
Eigen::MatrixXd matTObs
Definition sts_cluster.h:71
Eigen::MatrixXi matClusterIds
Definition sts_cluster.h:74
QVector< double > vecClusterStats
Definition sts_cluster.h:72
QVector< double > vecClusterPvals
Definition sts_cluster.h:73
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)