v2.0.0
Loading...
Searching...
No Matches
STSLIB::StatsCluster Class Reference

Cluster permutation testing. More...

#include <sts_cluster.h>

Static Public Member Functions

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 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 Eigen::MatrixXd tfce (const Eigen::MatrixXd &statMap, const Eigen::SparseMatrix< int > &adjacency, double E=0.5, double H=2.0, int nSteps=100)

Detailed Description

Cluster permutation testing.

Cluster-based permutation test for comparing two conditions.

Definition at line 84 of file sts_cluster.h.

Member Function Documentation

◆ fTestPermutationTest()

StatsClusterResult StatsCluster::fTestPermutationTest ( const QVector< QVector< Eigen::MatrixXd > > & conditions,
const Eigen::SparseMatrix< int > & adjacency,
double threshold,
int nPermutations )
static

F-test cluster-based permutation test for one-way ANOVA.

Tests for differences across conditions by randomly reassigning condition labels and computing max cluster F-statistics.

Parameters
[in]conditionsVector of conditions, each containing per-subject matrices (nVertices x nTimes).
[in]adjacencySpatio-temporal adjacency matrix (nVertices*nTimes x nVertices*nTimes).
[in]thresholdCluster-forming F-threshold.
[in]nPermutationsNumber of permutations.
Returns
StatsClusterResult with observed F-map (in matTObs), cluster statistics, and p-values.
Since
2.2.0

Definition at line 625 of file sts_cluster.cpp.

◆ oneSamplePermutationTest()

StatsClusterResult StatsCluster::oneSamplePermutationTest ( const QVector< Eigen::MatrixXd > & data,
const Eigen::SparseMatrix< int > & adjacency,
double threshold,
int nPermutations,
StatsTailType tail )
static

One-sample cluster-based permutation test.

Tests whether the mean across subjects differs from zero at each vertex x time point using sign-flip permutations.

Parameters
[in]dataPer-subject data. Each matrix is nVertices x nTimes.
[in]adjacencySpatio-temporal adjacency matrix (nVertices*nTimes x nVertices*nTimes).
[in]thresholdCluster-forming t-threshold.
[in]nPermutationsNumber of sign-flip permutations.
[in]tailTail type for the test.
Returns
StatsClusterResult with observed t-map, cluster statistics, and p-values.
Since
2.2.0

Definition at line 554 of file sts_cluster.cpp.

◆ permutationTest()

StatsClusterResult StatsCluster::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

Cluster-based permutation test.

Parameters
[in]dataAPer-subject data for condition A. Each matrix is nChannels x nTimes.
[in]dataBPer-subject data for condition B. Each matrix is nChannels x nTimes.
[in]adjacencySpatial adjacency matrix (nChannels x nChannels).
[in]nPermutationsNumber of permutations (default 1024).
[in]clusterAlphaAlpha level for initial thresholding (default 0.05).
[in]pThresholdp-value threshold for reporting significant clusters (default 0.05).
[in]tailTail type for the test.
Returns
StatsClusterResult with observed t-map, cluster statistics, and p-values.

Definition at line 62 of file sts_cluster.cpp.

◆ tfce()

MatrixXd StatsCluster::tfce ( const Eigen::MatrixXd & statMap,
const Eigen::SparseMatrix< int > & adjacency,
double E = 0.5,
double H = 2.0,
int nSteps = 100 )
static

Threshold-Free Cluster Enhancement (TFCE).

Enhances a statistic map by integrating cluster extent and height over a range of thresholds (Smith & Nichols 2009). Handles both positive and negative values.

Parameters
[in]statMapStatistic map (nVertices x nTimes).
[in]adjacencySpatio-temporal adjacency matrix (nVertices*nTimes x nVertices*nTimes).
[in]EExtent exponent (default 0.5).
[in]HHeight exponent (default 2.0).
[in]nStepsNumber of threshold steps (default 100).
Returns
TFCE-enhanced score map (nVertices x nTimes).
Since
2.2.0

Definition at line 688 of file sts_cluster.cpp.


The documentation for this class was generated from the following files: